Actual source code: dt.c
petsc-3.13.0 2020-03-29
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: const char *const PetscDTNodeTypes[] = {"gaussjacobi", "equispaced", "tanhsinh", "PETSCDTNODES_", 0};
17: static PetscBool GolubWelschCite = PETSC_FALSE;
18: const char GolubWelschCitation[] = "@article{GolubWelsch1969,\n"
19: " author = {Golub and Welsch},\n"
20: " title = {Calculation of Quadrature Rules},\n"
21: " journal = {Math. Comp.},\n"
22: " volume = {23},\n"
23: " number = {106},\n"
24: " pages = {221--230},\n"
25: " year = {1969}\n}\n";
27: /* Numerical tests in src/dm/dt/tests/ex1.c show that when computing the nodes and weights of Gauss-Jacobi
28: quadrature rules:
30: - in double precision, Newton's method and Golub & Welsch both work for moderate degrees (< 100),
31: - in single precision, Newton's method starts producing incorrect roots around n = 15, but
32: the weights from Golub & Welsch become a problem before then: they produces errors
33: in computing the Jacobi-polynomial Gram matrix around n = 6.
35: So we default to Newton's method (required fewer dependencies) */
36: PetscBool PetscDTGaussQuadratureNewton_Internal = PETSC_TRUE;
38: PetscClassId PETSCQUADRATURE_CLASSID = 0;
40: /*@
41: PetscQuadratureCreate - Create a PetscQuadrature object
43: Collective
45: Input Parameter:
46: . comm - The communicator for the PetscQuadrature object
48: Output Parameter:
49: . q - The PetscQuadrature object
51: Level: beginner
53: .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData()
54: @*/
55: PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
56: {
61: DMInitializePackage();
62: PetscHeaderCreate(*q,PETSCQUADRATURE_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);
63: (*q)->dim = -1;
64: (*q)->Nc = 1;
65: (*q)->order = -1;
66: (*q)->numPoints = 0;
67: (*q)->points = NULL;
68: (*q)->weights = NULL;
69: return(0);
70: }
72: /*@
73: PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object
75: Collective on q
77: Input Parameter:
78: . q - The PetscQuadrature object
80: Output Parameter:
81: . r - The new PetscQuadrature object
83: Level: beginner
85: .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData()
86: @*/
87: PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r)
88: {
89: PetscInt order, dim, Nc, Nq;
90: const PetscReal *points, *weights;
91: PetscReal *p, *w;
92: PetscErrorCode ierr;
96: PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);
97: PetscQuadratureGetOrder(q, &order);
98: PetscQuadratureSetOrder(*r, order);
99: PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);
100: PetscMalloc1(Nq*dim, &p);
101: PetscMalloc1(Nq*Nc, &w);
102: PetscArraycpy(p, points, Nq*dim);
103: PetscArraycpy(w, weights, Nc * Nq);
104: PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);
105: return(0);
106: }
108: /*@
109: PetscQuadratureDestroy - Destroys a PetscQuadrature object
111: Collective on q
113: Input Parameter:
114: . q - The PetscQuadrature object
116: Level: beginner
118: .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
119: @*/
120: PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
121: {
125: if (!*q) return(0);
127: if (--((PetscObject)(*q))->refct > 0) {
128: *q = NULL;
129: return(0);
130: }
131: PetscFree((*q)->points);
132: PetscFree((*q)->weights);
133: PetscHeaderDestroy(q);
134: return(0);
135: }
137: /*@
138: PetscQuadratureGetOrder - Return the order of the method
140: Not collective
142: Input Parameter:
143: . q - The PetscQuadrature object
145: Output Parameter:
146: . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
148: Level: intermediate
150: .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
151: @*/
152: PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order)
153: {
157: *order = q->order;
158: return(0);
159: }
161: /*@
162: PetscQuadratureSetOrder - Return the order of the method
164: Not collective
166: Input Parameters:
167: + q - The PetscQuadrature object
168: - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
170: Level: intermediate
172: .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
173: @*/
174: PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order)
175: {
178: q->order = order;
179: return(0);
180: }
182: /*@
183: PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated
185: Not collective
187: Input Parameter:
188: . q - The PetscQuadrature object
190: Output Parameter:
191: . Nc - The number of components
193: Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.
195: Level: intermediate
197: .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData()
198: @*/
199: PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc)
200: {
204: *Nc = q->Nc;
205: return(0);
206: }
208: /*@
209: PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated
211: Not collective
213: Input Parameters:
214: + q - The PetscQuadrature object
215: - Nc - The number of components
217: Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.
219: Level: intermediate
221: .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData()
222: @*/
223: PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc)
224: {
227: q->Nc = Nc;
228: return(0);
229: }
231: /*@C
232: PetscQuadratureGetData - Returns the data defining the quadrature
234: Not collective
236: Input Parameter:
237: . q - The PetscQuadrature object
239: Output Parameters:
240: + dim - The spatial dimension
241: . Nc - The number of components
242: . npoints - The number of quadrature points
243: . points - The coordinates of each quadrature point
244: - weights - The weight of each quadrature point
246: Level: intermediate
248: Fortran Notes:
249: From Fortran you must call PetscQuadratureRestoreData() when you are done with the data
251: .seealso: PetscQuadratureCreate(), PetscQuadratureSetData()
252: @*/
253: PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
254: {
257: if (dim) {
259: *dim = q->dim;
260: }
261: if (Nc) {
263: *Nc = q->Nc;
264: }
265: if (npoints) {
267: *npoints = q->numPoints;
268: }
269: if (points) {
271: *points = q->points;
272: }
273: if (weights) {
275: *weights = q->weights;
276: }
277: return(0);
278: }
280: static PetscErrorCode PetscDTJacobianInverse_Internal(PetscInt m, PetscInt n, const PetscReal J[], PetscReal Jinv[])
281: {
282: PetscScalar *Js, *Jinvs;
283: PetscInt i, j, k;
284: PetscBLASInt bm, bn, info;
288: if (!m || !n) return(0);
289: PetscBLASIntCast(m, &bm);
290: PetscBLASIntCast(n, &bn);
291: #if defined(PETSC_USE_COMPLEX)
292: PetscMalloc2(m*n, &Js, m*n, &Jinvs);
293: for (i = 0; i < m*n; i++) Js[i] = J[i];
294: #else
295: Js = (PetscReal *) J;
296: Jinvs = Jinv;
297: #endif
298: if (m == n) {
299: PetscBLASInt *pivots;
300: PetscScalar *W;
302: PetscMalloc2(m, &pivots, m, &W);
304: PetscArraycpy(Jinvs, Js, m * m);
305: PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, Jinvs, &bm, pivots, &info));
306: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info);
307: PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, Jinvs, &bm, pivots, W, &bm, &info));
308: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info);
309: PetscFree2(pivots, W);
310: } else if (m < n) {
311: PetscScalar *JJT;
312: PetscBLASInt *pivots;
313: PetscScalar *W;
315: PetscMalloc1(m*m, &JJT);
316: PetscMalloc2(m, &pivots, m, &W);
317: for (i = 0; i < m; i++) {
318: for (j = 0; j < m; j++) {
319: PetscScalar val = 0.;
321: for (k = 0; k < n; k++) val += Js[i * n + k] * Js[j * n + k];
322: JJT[i * m + j] = val;
323: }
324: }
326: PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, JJT, &bm, pivots, &info));
327: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info);
328: PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, JJT, &bm, pivots, W, &bm, &info));
329: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info);
330: for (i = 0; i < n; i++) {
331: for (j = 0; j < m; j++) {
332: PetscScalar val = 0.;
334: for (k = 0; k < m; k++) val += Js[k * n + i] * JJT[k * m + j];
335: Jinvs[i * m + j] = val;
336: }
337: }
338: PetscFree2(pivots, W);
339: PetscFree(JJT);
340: } else {
341: PetscScalar *JTJ;
342: PetscBLASInt *pivots;
343: PetscScalar *W;
345: PetscMalloc1(n*n, &JTJ);
346: PetscMalloc2(n, &pivots, n, &W);
347: for (i = 0; i < n; i++) {
348: for (j = 0; j < n; j++) {
349: PetscScalar val = 0.;
351: for (k = 0; k < m; k++) val += Js[k * n + i] * Js[k * n + j];
352: JTJ[i * n + j] = val;
353: }
354: }
356: PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bn, &bn, JTJ, &bn, pivots, &info));
357: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info);
358: PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bn, JTJ, &bn, pivots, W, &bn, &info));
359: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info);
360: for (i = 0; i < n; i++) {
361: for (j = 0; j < m; j++) {
362: PetscScalar val = 0.;
364: for (k = 0; k < n; k++) val += JTJ[i * n + k] * Js[j * n + k];
365: Jinvs[i * m + j] = val;
366: }
367: }
368: PetscFree2(pivots, W);
369: PetscFree(JTJ);
370: }
371: #if defined(PETSC_USE_COMPLEX)
372: for (i = 0; i < m*n; i++) Jinv[i] = PetscRealPart(Jinvs[i]);
373: PetscFree2(Js, Jinvs);
374: #endif
375: return(0);
376: }
378: /*@
379: PetscQuadraturePushForward - Push forward a quadrature functional under an affine transformation.
381: Collecive on PetscQuadrature
383: Input Arguments:
384: + q - the quadrature functional
385: . imageDim - the dimension of the image of the transformation
386: . origin - a point in the original space
387: . originImage - the image of the origin under the transformation
388: . J - the Jacobian of the image: an [imageDim x dim] matrix in row major order
389: - formDegree - transform the quadrature weights as k-forms of this form degree (if the number of components is a multiple of (dim choose formDegree), it is assumed that they represent multiple k-forms) [see PetscDTAltVPullback() for interpretation of formDegree]
391: Output Arguments:
392: . Jinvstarq - a quadrature rule where each point is the image of a point in the original quadrature rule, and where the k-form weights have been pulled-back by the pseudoinverse of J to the k-form weights in the image space.
394: Note: the new quadrature rule will have a different number of components if spaces have different dimensions. For example, pushing a 2-form forward from a two dimensional space to a three dimensional space changes the number of components from 1 to 3.
396: Level: intermediate
398: .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
399: @*/
400: PetscErrorCode PetscQuadraturePushForward(PetscQuadrature q, PetscInt imageDim, const PetscReal origin[], const PetscReal originImage[], const PetscReal J[], PetscInt formDegree, PetscQuadrature *Jinvstarq)
401: {
402: PetscInt dim, Nc, imageNc, formSize, Ncopies, imageFormSize, Npoints, pt, i, j, c;
403: const PetscReal *points;
404: const PetscReal *weights;
405: PetscReal *imagePoints, *imageWeights;
406: PetscReal *Jinv;
407: PetscReal *Jinvstar;
408: PetscErrorCode ierr;
412: if (imageDim < PetscAbsInt(formDegree)) SETERRQ2(PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Cannot represent a %D-form in %D dimensions", PetscAbsInt(formDegree), imageDim);
413: PetscQuadratureGetData(q, &dim, &Nc, &Npoints, &points, &weights);
414: PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &formSize);
415: if (Nc % formSize) SETERRQ2(PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Number of components %D is not a multiple of formSize %D\n", Nc, formSize);
416: Ncopies = Nc / formSize;
417: PetscDTBinomialInt(imageDim, PetscAbsInt(formDegree), &imageFormSize);
418: imageNc = Ncopies * imageFormSize;
419: PetscMalloc1(Npoints * imageDim, &imagePoints);
420: PetscMalloc1(Npoints * imageNc, &imageWeights);
421: PetscMalloc2(imageDim * dim, &Jinv, formSize * imageFormSize, &Jinvstar);
422: PetscDTJacobianInverse_Internal(imageDim, dim, J, Jinv);
423: PetscDTAltVPullbackMatrix(imageDim, dim, Jinv, formDegree, Jinvstar);
424: for (pt = 0; pt < Npoints; pt++) {
425: const PetscReal *point = &points[pt * dim];
426: PetscReal *imagePoint = &imagePoints[pt * imageDim];
428: for (i = 0; i < imageDim; i++) {
429: PetscReal val = originImage[i];
431: for (j = 0; j < dim; j++) val += J[i * dim + j] * (point[j] - origin[j]);
432: imagePoint[i] = val;
433: }
434: for (c = 0; c < Ncopies; c++) {
435: const PetscReal *form = &weights[pt * Nc + c * formSize];
436: PetscReal *imageForm = &imageWeights[pt * imageNc + c * imageFormSize];
438: for (i = 0; i < imageFormSize; i++) {
439: PetscReal val = 0.;
441: for (j = 0; j < formSize; j++) val += Jinvstar[i * formSize + j] * form[j];
442: imageForm[i] = val;
443: }
444: }
445: }
446: PetscQuadratureCreate(PetscObjectComm((PetscObject)q), Jinvstarq);
447: PetscQuadratureSetData(*Jinvstarq, imageDim, imageNc, Npoints, imagePoints, imageWeights);
448: PetscFree2(Jinv, Jinvstar);
449: return(0);
450: }
452: /*@C
453: PetscQuadratureSetData - Sets the data defining the quadrature
455: Not collective
457: Input Parameters:
458: + q - The PetscQuadrature object
459: . dim - The spatial dimension
460: . Nc - The number of components
461: . npoints - The number of quadrature points
462: . points - The coordinates of each quadrature point
463: - weights - The weight of each quadrature point
465: Note: This routine owns the references to points and weights, so they must be allocated using PetscMalloc() and the user should not free them.
467: Level: intermediate
469: .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
470: @*/
471: PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
472: {
475: if (dim >= 0) q->dim = dim;
476: if (Nc >= 0) q->Nc = Nc;
477: if (npoints >= 0) q->numPoints = npoints;
478: if (points) {
480: q->points = points;
481: }
482: if (weights) {
484: q->weights = weights;
485: }
486: return(0);
487: }
489: static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v)
490: {
491: PetscInt q, d, c;
492: PetscViewerFormat format;
493: PetscErrorCode ierr;
496: if (quad->Nc > 1) {PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D) with %D components\n", quad->order, quad->numPoints, quad->dim, quad->Nc);}
497: else {PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D)\n", quad->order, quad->numPoints, quad->dim);}
498: PetscViewerGetFormat(v, &format);
499: if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
500: for (q = 0; q < quad->numPoints; ++q) {
501: PetscViewerASCIIPrintf(v, "p%D (", q);
502: PetscViewerASCIIUseTabs(v, PETSC_FALSE);
503: for (d = 0; d < quad->dim; ++d) {
504: if (d) {PetscViewerASCIIPrintf(v, ", ");}
505: PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q*quad->dim+d]);
506: }
507: PetscViewerASCIIPrintf(v, ") ");
508: if (quad->Nc > 1) {PetscViewerASCIIPrintf(v, "w%D (", q);}
509: for (c = 0; c < quad->Nc; ++c) {
510: if (c) {PetscViewerASCIIPrintf(v, ", ");}
511: PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q*quad->Nc+c]);
512: }
513: if (quad->Nc > 1) {PetscViewerASCIIPrintf(v, ")");}
514: PetscViewerASCIIPrintf(v, "\n");
515: PetscViewerASCIIUseTabs(v, PETSC_TRUE);
516: }
517: return(0);
518: }
520: /*@C
521: PetscQuadratureView - Views a PetscQuadrature object
523: Collective on quad
525: Input Parameters:
526: + quad - The PetscQuadrature object
527: - viewer - The PetscViewer object
529: Level: beginner
531: .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
532: @*/
533: PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
534: {
535: PetscBool iascii;
541: if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) quad), &viewer);}
542: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
543: PetscViewerASCIIPushTab(viewer);
544: if (iascii) {PetscQuadratureView_Ascii(quad, viewer);}
545: PetscViewerASCIIPopTab(viewer);
546: return(0);
547: }
549: /*@C
550: PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement
552: Not collective
554: Input Parameter:
555: + q - The original PetscQuadrature
556: . numSubelements - The number of subelements the original element is divided into
557: . v0 - An array of the initial points for each subelement
558: - jac - An array of the Jacobian mappings from the reference to each subelement
560: Output Parameters:
561: . dim - The dimension
563: Note: Together v0 and jac define an affine mapping from the original reference element to each subelement
565: Not available from Fortran
567: Level: intermediate
569: .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
570: @*/
571: PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref)
572: {
573: const PetscReal *points, *weights;
574: PetscReal *pointsRef, *weightsRef;
575: PetscInt dim, Nc, order, npoints, npointsRef, c, p, cp, d, e;
576: PetscErrorCode ierr;
583: PetscQuadratureCreate(PETSC_COMM_SELF, qref);
584: PetscQuadratureGetOrder(q, &order);
585: PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);
586: npointsRef = npoints*numSubelements;
587: PetscMalloc1(npointsRef*dim,&pointsRef);
588: PetscMalloc1(npointsRef*Nc, &weightsRef);
589: for (c = 0; c < numSubelements; ++c) {
590: for (p = 0; p < npoints; ++p) {
591: for (d = 0; d < dim; ++d) {
592: pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d];
593: for (e = 0; e < dim; ++e) {
594: pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0);
595: }
596: }
597: /* Could also use detJ here */
598: for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements;
599: }
600: }
601: PetscQuadratureSetOrder(*qref, order);
602: PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);
603: return(0);
604: }
606: /* Compute the coefficients for the Jacobi polynomial recurrence,
607: *
608: * J^{a,b}_n(x) = (cnm1 + cnm1x * x) * J^{a,b}_{n-1}(x) - cnm2 * J^{a,b}_{n-2}(x).
609: */
610: #define PetscDTJacobiRecurrence_Internal(n,a,b,cnm1,cnm1x,cnm2) \
611: do { \
612: PetscReal _a = (a); \
613: PetscReal _b = (b); \
614: PetscReal _n = (n); \
615: if (n == 1) { \
616: (cnm1) = (_a-_b) * 0.5; \
617: (cnm1x) = (_a+_b+2.)*0.5; \
618: (cnm2) = 0.; \
619: } else { \
620: PetscReal _2n = _n+_n; \
621: PetscReal _d = (_2n*(_n+_a+_b)*(_2n+_a+_b-2)); \
622: PetscReal _n1 = (_2n+_a+_b-1.)*(_a*_a-_b*_b); \
623: PetscReal _n1x = (_2n+_a+_b-1.)*(_2n+_a+_b)*(_2n+_a+_b-2); \
624: PetscReal _n2 = 2.*((_n+_a-1.)*(_n+_b-1.)*(_2n+_a+_b)); \
625: (cnm1) = _n1 / _d; \
626: (cnm1x) = _n1x / _d; \
627: (cnm2) = _n2 / _d; \
628: } \
629: } while (0)
631: static PetscErrorCode PetscDTJacobiEval_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscInt k, const PetscReal *points, PetscInt ndegree, const PetscInt *degrees, PetscReal *p)
632: {
633: PetscReal ak, bk;
634: PetscReal abk1;
635: PetscInt i,l,maxdegree;
638: maxdegree = degrees[ndegree-1] - k;
639: ak = a + k;
640: bk = b + k;
641: abk1 = a + b + k + 1.;
642: if (maxdegree < 0) {
643: for (i = 0; i < npoints; i++) for (l = 0; l < ndegree; l++) p[i*ndegree+l] = 0.;
644: return(0);
645: }
646: for (i=0; i<npoints; i++) {
647: PetscReal pm1,pm2,x;
648: PetscReal cnm1, cnm1x, cnm2;
649: PetscInt j,m;
651: x = points[i];
652: pm2 = 1.;
653: PetscDTJacobiRecurrence_Internal(1,ak,bk,cnm1,cnm1x,cnm2);
654: pm1 = (cnm1 + cnm1x*x);
655: l = 0;
656: while (l < ndegree && degrees[l] - k < 0) {
657: p[l++] = 0.;
658: }
659: while (l < ndegree && degrees[l] - k == 0) {
660: p[l] = pm2;
661: for (m = 0; m < k; m++) p[l] *= (abk1 + m) * 0.5;
662: l++;
663: }
664: while (l < ndegree && degrees[l] - k == 1) {
665: p[l] = pm1;
666: for (m = 0; m < k; m++) p[l] *= (abk1 + 1 + m) * 0.5;
667: l++;
668: }
669: for (j=2; j<=maxdegree; j++) {
670: PetscReal pp;
672: PetscDTJacobiRecurrence_Internal(j,ak,bk,cnm1,cnm1x,cnm2);
673: pp = (cnm1 + cnm1x*x)*pm1 - cnm2*pm2;
674: pm2 = pm1;
675: pm1 = pp;
676: while (l < ndegree && degrees[l] - k == j) {
677: p[l] = pp;
678: for (m = 0; m < k; m++) p[l] *= (abk1 + j + m) * 0.5;
679: l++;
680: }
681: }
682: p += ndegree;
683: }
684: return(0);
685: }
687: /*@
688: PetscDTJacobiEval - evaluate Jacobi polynomials for the weight function $(1.+x)^{\alpha} (1.-x)^{\beta}$
689: at points
691: Not Collective
693: Input Arguments:
694: + npoints - number of spatial points to evaluate at
695: . alpha - the left exponent > -1
696: . beta - the right exponent > -1
697: . points - array of locations to evaluate at
698: . ndegree - number of basis degrees to evaluate
699: - degrees - sorted array of degrees to evaluate
701: Output Arguments:
702: + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
703: . D - row-oriented derivative evaluation matrix (or NULL)
704: - D2 - row-oriented second derivative evaluation matrix (or NULL)
706: Level: intermediate
708: .seealso: PetscDTGaussQuadrature()
709: @*/
710: PetscErrorCode PetscDTJacobiEval(PetscInt npoints,PetscReal alpha, PetscReal beta, const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
711: {
715: if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1.");
716: if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1.");
717: if (!npoints || !ndegree) return(0);
718: if (B) {PetscDTJacobiEval_Internal(npoints, alpha, beta, 0, points, ndegree, degrees, B);}
719: if (D) {PetscDTJacobiEval_Internal(npoints, alpha, beta, 1, points, ndegree, degrees, D);}
720: if (D2) {PetscDTJacobiEval_Internal(npoints, alpha, beta, 2, points, ndegree, degrees, D2);}
721: return(0);
722: }
724: /*@
725: PetscDTLegendreEval - evaluate Legendre polynomials at points
727: Not Collective
729: Input Arguments:
730: + npoints - number of spatial points to evaluate at
731: . points - array of locations to evaluate at
732: . ndegree - number of basis degrees to evaluate
733: - degrees - sorted array of degrees to evaluate
735: Output Arguments:
736: + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
737: . D - row-oriented derivative evaluation matrix (or NULL)
738: - D2 - row-oriented second derivative evaluation matrix (or NULL)
740: Level: intermediate
742: .seealso: PetscDTGaussQuadrature()
743: @*/
744: PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
745: {
749: PetscDTJacobiEval(npoints, 0., 0., points, ndegree, degrees, B, D, D2);
750: return(0);
751: }
753: /* solve the symmetric tridiagonal eigenvalue system, writing the eigenvalues into eigs and the eigenvectors into V
754: * with lds n; diag and subdiag are overwritten */
755: static PetscErrorCode PetscDTSymmetricTridiagonalEigensolve(PetscInt n, PetscReal diag[], PetscReal subdiag[],
756: PetscReal eigs[], PetscScalar V[])
757: {
758: char jobz = 'V'; /* eigenvalues and eigenvectors */
759: char range = 'A'; /* all eigenvalues will be found */
760: PetscReal VL = 0.; /* ignored because range is 'A' */
761: PetscReal VU = 0.; /* ignored because range is 'A' */
762: PetscBLASInt IL = 0; /* ignored because range is 'A' */
763: PetscBLASInt IU = 0; /* ignored because range is 'A' */
764: PetscReal abstol = 0.; /* unused */
765: PetscBLASInt bn, bm, ldz; /* bm will equal bn on exit */
766: PetscBLASInt *isuppz;
767: PetscBLASInt lwork, liwork;
768: PetscReal workquery;
769: PetscBLASInt iworkquery;
770: PetscBLASInt *iwork;
771: PetscBLASInt info;
772: PetscReal *work = NULL;
776: #if !defined(PETSCDTGAUSSIANQUADRATURE_EIG)
777: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found");
778: #endif
779: PetscBLASIntCast(n, &bn);
780: PetscBLASIntCast(n, &ldz);
781: #if !defined(PETSC_MISSING_LAPACK_STEGR)
782: PetscMalloc1(2 * n, &isuppz);
783: lwork = -1;
784: liwork = -1;
785: PetscStackCallBLAS("LAPACKstegr",LAPACKstegr_(&jobz,&range,&bn,diag,subdiag,&VL,&VU,&IL,&IU,&abstol,&bm,eigs,V,&ldz,isuppz,&workquery,&lwork,&iworkquery,&liwork,&info));
786: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEGR error");
787: lwork = (PetscBLASInt) workquery;
788: liwork = (PetscBLASInt) iworkquery;
789: PetscMalloc2(lwork, &work, liwork, &iwork);
790: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
791: PetscStackCallBLAS("LAPACKstegr",LAPACKstegr_(&jobz,&range,&bn,diag,subdiag,&VL,&VU,&IL,&IU,&abstol,&bm,eigs,V,&ldz,isuppz,work,&lwork,iwork,&liwork,&info));
792: PetscFPTrapPop();
793: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEGR error");
794: PetscFree2(work, iwork);
795: PetscFree(isuppz);
796: #elif !defined(PETSC_MISSING_LAPACK_STEQR)
797: jobz = 'I'; /* Compute eigenvalues and eigenvectors of the
798: tridiagonal matrix. Z is initialized to the identity
799: matrix. */
800: PetscMalloc1(PetscMax(1,2*n-2),&work);
801: PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&bn,diag,subdiag,V,&ldz,work,&info));
802: PetscFPTrapPop();
803: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
804: PetscFree(work);
805: PetscArraycpy(eigs,diag,n);
806: #endif
807: return(0);
808: }
810: /* Formula for the weights at the endpoints (-1 and 1) of Gauss-Lobatto-Jacobi
811: * quadrature rules on the interval [-1, 1] */
812: static PetscErrorCode PetscDTGaussLobattoJacobiEndweights_Internal(PetscInt n, PetscReal alpha, PetscReal beta, PetscReal *leftw, PetscReal *rightw)
813: {
814: PetscReal twoab1;
815: PetscInt m = n - 2;
816: PetscReal a = alpha + 1.;
817: PetscReal b = beta + 1.;
818: PetscReal gra, grb;
821: twoab1 = PetscPowReal(2., a + b - 1.);
822: #if defined(PETSC_HAVE_LGAMMA)
823: grb = PetscExpReal(2. * PetscLGamma(b+1.) + PetscLGamma(m+1.) + PetscLGamma(m+a+1.) -
824: (PetscLGamma(m+b+1) + PetscLGamma(m+a+b+1.)));
825: gra = PetscExpReal(2. * PetscLGamma(a+1.) + PetscLGamma(m+1.) + PetscLGamma(m+b+1.) -
826: (PetscLGamma(m+a+1) + PetscLGamma(m+a+b+1.)));
827: #else
828: {
829: PetscInt alphai = (PetscInt) alpha;
830: PetscInt betai = (PetscInt) beta;
833: if ((PetscReal) alphai == alpha && (PetscReal) betai == beta) {
834: PetscReal binom1, binom2;
836: PetscDTBinomial(m+b, b, &binom1);
837: PetscDTBinomial(m+a+b, b, &binom2);
838: grb = 1./ (binom1 * binom2);
839: PetscDTBinomial(m+a, a, &binom1);
840: PetscDTBinomial(m+a+b, a, &binom2);
841: gra = 1./ (binom1 * binom2);
842: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable.");
843: }
844: #endif
845: *leftw = twoab1 * grb / b;
846: *rightw = twoab1 * gra / a;
847: return(0);
848: }
850: /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
851: Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
852: PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
853: {
854: PetscReal pn1, pn2;
855: PetscReal cnm1, cnm1x, cnm2;
856: PetscInt k;
859: if (!n) {*P = 1.0; return(0);}
860: PetscDTJacobiRecurrence_Internal(1,a,b,cnm1,cnm1x,cnm2);
861: pn2 = 1.;
862: pn1 = cnm1 + cnm1x*x;
863: if (n == 1) {*P = pn1; return(0);}
864: *P = 0.0;
865: for (k = 2; k < n+1; ++k) {
866: PetscDTJacobiRecurrence_Internal(k,a,b,cnm1,cnm1x,cnm2);
868: *P = (cnm1 + cnm1x*x)*pn1 - cnm2*pn2;
869: pn2 = pn1;
870: pn1 = *P;
871: }
872: return(0);
873: }
875: /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
876: PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscInt k, PetscReal *P)
877: {
878: PetscReal nP;
879: PetscInt i;
883: if (k > n) {*P = 0.0; return(0);}
884: PetscDTComputeJacobi(a+k, b+k, n-k, x, &nP);
885: for (i = 0; i < k; i++) nP *= (a + b + n + 1. + i) * 0.5;
886: *P = nP;
887: return(0);
888: }
890: /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
891: PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
892: {
894: *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
895: *eta = y;
896: return(0);
897: }
899: static PetscErrorCode PetscDTGaussJacobiQuadrature_Newton_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal x[], PetscReal w[])
900: {
901: PetscInt maxIter = 100;
902: PetscReal eps = PetscExpReal(0.75 * PetscLogReal(PETSC_MACHINE_EPSILON));
903: PetscReal a1, a6, gf;
904: PetscInt k;
909: a1 = PetscPowReal(2.0, a+b+1);
910: #if defined(PETSC_HAVE_LGAMMA)
911: {
912: PetscReal a2, a3, a4, a5;
913: a2 = PetscLGamma(a + npoints + 1);
914: a3 = PetscLGamma(b + npoints + 1);
915: a4 = PetscLGamma(a + b + npoints + 1);
916: a5 = PetscLGamma(npoints + 1);
917: gf = PetscExpReal(a2 + a3 - (a4 + a5));
918: }
919: #else
920: {
921: PetscInt ia, ib;
923: ia = (PetscInt) a;
924: ib = (PetscInt) b;
925: gf = 1.;
926: if (ia == a && ia >= 0) { /* compute ratio of rising factorals wrt a */
927: for (k = 0; k < ia; k++) gf *= (npoints + 1. + k) / (npoints + b + 1. + k);
928: } else if (b == b && ib >= 0) { /* compute ratio of rising factorials wrt b */
929: for (k = 0; k < ib; k++) gf *= (npoints + 1. + k) / (npoints + a + 1. + k);
930: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable.");
931: }
932: #endif
934: a6 = a1 * gf;
935: /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
936: Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
937: for (k = 0; k < npoints; ++k) {
938: PetscReal r = PetscCosReal(PETSC_PI * (1. - (4.*k + 3. + 2.*b) / (4.*npoints + 2.*(a + b + 1.)))), dP;
939: PetscInt j;
941: if (k > 0) r = 0.5 * (r + x[k-1]);
942: for (j = 0; j < maxIter; ++j) {
943: PetscReal s = 0.0, delta, f, fp;
944: PetscInt i;
946: for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
947: PetscDTComputeJacobi(a, b, npoints, r, &f);
948: PetscDTComputeJacobiDerivative(a, b, npoints, r, 1, &fp);
949: delta = f / (fp - f * s);
950: r = r - delta;
951: if (PetscAbsReal(delta) < eps) break;
952: }
953: x[k] = r;
954: PetscDTComputeJacobiDerivative(a, b, npoints, x[k], 1, &dP);
955: w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
956: }
957: return(0);
958: }
960: /* Compute the diagonals of the Jacobi matrix used in Golub & Welsch algorithms for Gauss-Jacobi
961: * quadrature weight calculations on [-1,1] for exponents (1. + x)^a (1.-x)^b */
962: static PetscErrorCode PetscDTJacobiMatrix_Internal(PetscInt nPoints, PetscReal a, PetscReal b, PetscReal *d, PetscReal *s)
963: {
964: PetscInt i;
967: for (i = 0; i < nPoints; i++) {
968: PetscReal A, B, C;
970: PetscDTJacobiRecurrence_Internal(i+1,a,b,A,B,C);
971: d[i] = -A / B;
972: if (i) s[i-1] *= C / B;
973: if (i < nPoints - 1) s[i] = 1. / B;
974: }
975: return(0);
976: }
978: static PetscErrorCode PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
979: {
980: PetscReal mu0;
981: PetscReal ga, gb, gab;
982: PetscInt i;
986: PetscCitationsRegister(GolubWelschCitation, &GolubWelschCite);
988: #if defined(PETSC_HAVE_TGAMMA)
989: ga = PetscTGamma(a + 1);
990: gb = PetscTGamma(b + 1);
991: gab = PetscTGamma(a + b + 2);
992: #else
993: {
994: PetscInt ia, ib;
996: ia = (PetscInt) a;
997: ib = (PetscInt) b;
998: if (ia == a && ib == b && ia + 1 > 0 && ib + 1 > 0 && ia + ib + 2 > 0) { /* All gamma(x) terms are (x-1)! terms */
999: PetscDTFactorial(ia, &ga);
1000: PetscDTFactorial(ib, &gb);
1001: PetscDTFactorial(ia + ib + 1, &gb);
1002: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
1003: }
1004: #endif
1005: mu0 = PetscPowReal(2.,a + b + 1.) * ga * gb / gab;
1007: #if defined(PETSCDTGAUSSIANQUADRATURE_EIG)
1008: {
1009: PetscReal *diag, *subdiag;
1010: PetscScalar *V;
1012: PetscMalloc2(npoints, &diag, npoints, &subdiag);
1013: PetscMalloc1(npoints*npoints, &V);
1014: PetscDTJacobiMatrix_Internal(npoints, a, b, diag, subdiag);
1015: for (i = 0; i < npoints - 1; i++) subdiag[i] = PetscSqrtReal(subdiag[i]);
1016: PetscDTSymmetricTridiagonalEigensolve(npoints, diag, subdiag, x, V);
1017: for (i = 0; i < npoints; i++) w[i] = PetscSqr(PetscRealPart(V[i * npoints])) * mu0;
1018: PetscFree(V);
1019: PetscFree2(diag, subdiag);
1020: }
1021: #else
1022: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found");
1023: #endif
1024: { /* As of March 2, 2020, The Sun Performance Library breaks the LAPACK contract for xstegr and xsteqr: the
1025: eigenvalues are not guaranteed to be in ascending order. So we heave a passive aggressive sigh and check that
1026: the eigenvalues are sorted */
1027: PetscBool sorted;
1029: PetscSortedReal(npoints, x, &sorted);
1030: if (!sorted) {
1031: PetscInt *order, i;
1032: PetscReal *tmp;
1034: PetscMalloc2(npoints, &order, npoints, &tmp);
1035: for (i = 0; i < npoints; i++) order[i] = i;
1036: PetscSortRealWithPermutation(npoints, x, order);
1037: PetscArraycpy(tmp, x, npoints);
1038: for (i = 0; i < npoints; i++) x[i] = tmp[order[i]];
1039: PetscArraycpy(tmp, w, npoints);
1040: for (i = 0; i < npoints; i++) w[i] = tmp[order[i]];
1041: PetscFree2(order, tmp);
1042: }
1043: }
1044: return(0);
1045: }
1047: static PetscErrorCode PetscDTGaussJacobiQuadrature_Internal(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton)
1048: {
1052: if (npoints < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be positive");
1053: /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */
1054: if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1.");
1055: if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1.");
1057: if (newton) {
1058: PetscDTGaussJacobiQuadrature_Newton_Internal(npoints, alpha, beta, x, w);
1059: } else {
1060: PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(npoints, alpha, beta, x, w);
1061: }
1062: if (alpha == beta) { /* symmetrize */
1063: PetscInt i;
1064: for (i = 0; i < (npoints + 1) / 2; i++) {
1065: PetscInt j = npoints - 1 - i;
1066: PetscReal xi = x[i];
1067: PetscReal xj = x[j];
1068: PetscReal wi = w[i];
1069: PetscReal wj = w[j];
1071: x[i] = (xi - xj) / 2.;
1072: x[j] = (xj - xi) / 2.;
1073: w[i] = w[j] = (wi + wj) / 2.;
1074: }
1075: }
1076: return(0);
1077: }
1079: /*@
1080: PetscDTGaussJacobiQuadrature - quadrature for the interval [a, b] with the weight function
1081: $(x-a)^\alpha (x-b)^\beta$.
1083: Not collective
1085: Input Parameters:
1086: + npoints - the number of points in the quadrature rule
1087: . a - the left endpoint of the interval
1088: . b - the right endpoint of the interval
1089: . alpha - the left exponent
1090: - beta - the right exponent
1092: Output Parameters:
1093: + x - array of length npoints, the locations of the quadrature points
1094: - w - array of length npoints, the weights of the quadrature points
1096: Level: intermediate
1098: Note: this quadrature rule is exact for polynomials up to degree 2*npoints - 1.
1099: @*/
1100: PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt npoints,PetscReal a, PetscReal b, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[])
1101: {
1102: PetscInt i;
1106: PetscDTGaussJacobiQuadrature_Internal(npoints, alpha, beta, x, w, PetscDTGaussQuadratureNewton_Internal);
1107: if (a != -1. || b != 1.) { /* shift */
1108: for (i = 0; i < npoints; i++) {
1109: x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
1110: w[i] *= (b - a) / 2.;
1111: }
1112: }
1113: return(0);
1114: }
1116: static PetscErrorCode PetscDTGaussLobattoJacobiQuadrature_Internal(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton)
1117: {
1118: PetscInt i;
1122: if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be positive");
1123: /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */
1124: if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1.");
1125: if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1.");
1127: x[0] = -1.;
1128: x[npoints-1] = 1.;
1129: if (npoints > 2) {
1130: PetscDTGaussJacobiQuadrature_Internal(npoints-2, alpha+1., beta+1., &x[1], &w[1], newton);
1131: }
1132: for (i = 1; i < npoints - 1; i++) {
1133: w[i] /= (1. - x[i]*x[i]);
1134: }
1135: PetscDTGaussLobattoJacobiEndweights_Internal(npoints, alpha, beta, &w[0], &w[npoints-1]);
1136: return(0);
1137: }
1139: /*@
1140: PetscDTGaussLobattoJacobiQuadrature - quadrature for the interval [a, b] with the weight function
1141: $(x-a)^\alpha (x-b)^\beta$, with endpoints a and b included as quadrature points.
1143: Not collective
1145: Input Parameters:
1146: + npoints - the number of points in the quadrature rule
1147: . a - the left endpoint of the interval
1148: . b - the right endpoint of the interval
1149: . alpha - the left exponent
1150: - beta - the right exponent
1152: Output Parameters:
1153: + x - array of length npoints, the locations of the quadrature points
1154: - w - array of length npoints, the weights of the quadrature points
1156: Level: intermediate
1158: Note: this quadrature rule is exact for polynomials up to degree 2*npoints - 3.
1159: @*/
1160: PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt npoints,PetscReal a, PetscReal b, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[])
1161: {
1162: PetscInt i;
1166: PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, alpha, beta, x, w, PetscDTGaussQuadratureNewton_Internal);
1167: if (a != -1. || b != 1.) { /* shift */
1168: for (i = 0; i < npoints; i++) {
1169: x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
1170: w[i] *= (b - a) / 2.;
1171: }
1172: }
1173: return(0);
1174: }
1176: /*@
1177: PetscDTGaussQuadrature - create Gauss-Legendre quadrature
1179: Not Collective
1181: Input Arguments:
1182: + npoints - number of points
1183: . a - left end of interval (often-1)
1184: - b - right end of interval (often +1)
1186: Output Arguments:
1187: + x - quadrature points
1188: - w - quadrature weights
1190: Level: intermediate
1192: References:
1193: . 1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969.
1195: .seealso: PetscDTLegendreEval()
1196: @*/
1197: PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
1198: {
1199: PetscInt i;
1203: PetscDTGaussJacobiQuadrature_Internal(npoints, 0., 0., x, w, PetscDTGaussQuadratureNewton_Internal);
1204: if (a != -1. || b != 1.) { /* shift */
1205: for (i = 0; i < npoints; i++) {
1206: x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
1207: w[i] *= (b - a) / 2.;
1208: }
1209: }
1210: return(0);
1211: }
1213: /*@C
1214: PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre
1215: nodes of a given size on the domain [-1,1]
1217: Not Collective
1219: Input Parameter:
1220: + n - number of grid nodes
1221: - type - PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA or PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON
1223: Output Arguments:
1224: + x - quadrature points
1225: - w - quadrature weights
1227: Notes:
1228: For n > 30 the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not
1229: close enough to the desired solution
1231: These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes
1233: See https://epubs.siam.org/doi/abs/10.1137/110855442 https://epubs.siam.org/doi/abs/10.1137/120889873 for better ways to compute GLL nodes
1235: Level: intermediate
1237: .seealso: PetscDTGaussQuadrature()
1239: @*/
1240: PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w)
1241: {
1242: PetscBool newton;
1246: if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must provide at least 2 grid points per element");
1247: newton = (PetscBool) (type == PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON);
1248: PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, 0., 0., x, w, newton);
1249: return(0);
1250: }
1252: /*@
1253: PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature
1255: Not Collective
1257: Input Arguments:
1258: + dim - The spatial dimension
1259: . Nc - The number of components
1260: . npoints - number of points in one dimension
1261: . a - left end of interval (often-1)
1262: - b - right end of interval (often +1)
1264: Output Argument:
1265: . q - A PetscQuadrature object
1267: Level: intermediate
1269: .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval()
1270: @*/
1271: PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
1272: {
1273: PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c;
1274: PetscReal *x, *w, *xw, *ww;
1278: PetscMalloc1(totpoints*dim,&x);
1279: PetscMalloc1(totpoints*Nc,&w);
1280: /* Set up the Golub-Welsch system */
1281: switch (dim) {
1282: case 0:
1283: PetscFree(x);
1284: PetscFree(w);
1285: PetscMalloc1(1, &x);
1286: PetscMalloc1(Nc, &w);
1287: x[0] = 0.0;
1288: for (c = 0; c < Nc; ++c) w[c] = 1.0;
1289: break;
1290: case 1:
1291: PetscMalloc1(npoints,&ww);
1292: PetscDTGaussQuadrature(npoints, a, b, x, ww);
1293: for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i];
1294: PetscFree(ww);
1295: break;
1296: case 2:
1297: PetscMalloc2(npoints,&xw,npoints,&ww);
1298: PetscDTGaussQuadrature(npoints, a, b, xw, ww);
1299: for (i = 0; i < npoints; ++i) {
1300: for (j = 0; j < npoints; ++j) {
1301: x[(i*npoints+j)*dim+0] = xw[i];
1302: x[(i*npoints+j)*dim+1] = xw[j];
1303: for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j];
1304: }
1305: }
1306: PetscFree2(xw,ww);
1307: break;
1308: case 3:
1309: PetscMalloc2(npoints,&xw,npoints,&ww);
1310: PetscDTGaussQuadrature(npoints, a, b, xw, ww);
1311: for (i = 0; i < npoints; ++i) {
1312: for (j = 0; j < npoints; ++j) {
1313: for (k = 0; k < npoints; ++k) {
1314: x[((i*npoints+j)*npoints+k)*dim+0] = xw[i];
1315: x[((i*npoints+j)*npoints+k)*dim+1] = xw[j];
1316: x[((i*npoints+j)*npoints+k)*dim+2] = xw[k];
1317: for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k];
1318: }
1319: }
1320: }
1321: PetscFree2(xw,ww);
1322: break;
1323: default:
1324: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
1325: }
1326: PetscQuadratureCreate(PETSC_COMM_SELF, q);
1327: PetscQuadratureSetOrder(*q, 2*npoints-1);
1328: PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);
1329: PetscObjectChangeTypeName((PetscObject)*q,"GaussTensor");
1330: return(0);
1331: }
1333: /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
1334: PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
1335: {
1337: *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
1338: *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0;
1339: *zeta = z;
1340: return(0);
1341: }
1344: /*@
1345: PetscDTStroudConicalQuadrature - create Stroud conical quadrature for a simplex
1347: Not Collective
1349: Input Arguments:
1350: + dim - The simplex dimension
1351: . Nc - The number of components
1352: . npoints - The number of points in one dimension
1353: . a - left end of interval (often-1)
1354: - b - right end of interval (often +1)
1356: Output Argument:
1357: . q - A PetscQuadrature object
1359: Level: intermediate
1361: References:
1362: . 1. - Karniadakis and Sherwin. FIAT
1364: Note: For dim == 1, this is Gauss-Legendre quadrature
1366: .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature()
1367: @*/
1368: PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
1369: {
1370: PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints;
1371: PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w;
1375: if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
1376: PetscMalloc1(totpoints*dim, &x);
1377: PetscMalloc1(totpoints*Nc, &w);
1378: switch (dim) {
1379: case 0:
1380: PetscFree(x);
1381: PetscFree(w);
1382: PetscMalloc1(1, &x);
1383: PetscMalloc1(Nc, &w);
1384: x[0] = 0.0;
1385: for (c = 0; c < Nc; ++c) w[c] = 1.0;
1386: break;
1387: case 1:
1388: PetscMalloc1(npoints,&wx);
1389: PetscDTGaussJacobiQuadrature(npoints, -1., 1., 0.0, 0.0, x, wx);
1390: for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i];
1391: PetscFree(wx);
1392: break;
1393: case 2:
1394: PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);
1395: PetscDTGaussJacobiQuadrature(npoints, -1., 1., 0.0, 0.0, px, wx);
1396: PetscDTGaussJacobiQuadrature(npoints, -1., 1., 1.0, 0.0, py, wy);
1397: for (i = 0; i < npoints; ++i) {
1398: for (j = 0; j < npoints; ++j) {
1399: PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);
1400: for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j];
1401: }
1402: }
1403: PetscFree4(px,wx,py,wy);
1404: break;
1405: case 3:
1406: PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);
1407: PetscDTGaussJacobiQuadrature(npoints, -1., 1., 0.0, 0.0, px, wx);
1408: PetscDTGaussJacobiQuadrature(npoints, -1., 1., 1.0, 0.0, py, wy);
1409: PetscDTGaussJacobiQuadrature(npoints, -1., 1., 2.0, 0.0, pz, wz);
1410: for (i = 0; i < npoints; ++i) {
1411: for (j = 0; j < npoints; ++j) {
1412: for (k = 0; k < npoints; ++k) {
1413: 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]);
1414: for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k];
1415: }
1416: }
1417: }
1418: PetscFree6(px,wx,py,wy,pz,wz);
1419: break;
1420: default:
1421: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
1422: }
1423: PetscQuadratureCreate(PETSC_COMM_SELF, q);
1424: PetscQuadratureSetOrder(*q, 2*npoints-1);
1425: PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);
1426: PetscObjectChangeTypeName((PetscObject)*q,"GaussJacobi");
1427: return(0);
1428: }
1430: /*@
1431: PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell
1433: Not Collective
1435: Input Arguments:
1436: + dim - The cell dimension
1437: . level - The number of points in one dimension, 2^l
1438: . a - left end of interval (often-1)
1439: - b - right end of interval (often +1)
1441: Output Argument:
1442: . q - A PetscQuadrature object
1444: Level: intermediate
1446: .seealso: PetscDTGaussTensorQuadrature()
1447: @*/
1448: PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q)
1449: {
1450: const PetscInt p = 16; /* Digits of precision in the evaluation */
1451: const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */
1452: const PetscReal beta = (b+a)/2.; /* Center of the integration interval */
1453: const PetscReal h = PetscPowReal(2.0, -level); /* Step size, length between x_k */
1454: PetscReal xk; /* Quadrature point x_k on reference domain [-1, 1] */
1455: PetscReal wk = 0.5*PETSC_PI; /* Quadrature weight at x_k */
1456: PetscReal *x, *w;
1457: PetscInt K, k, npoints;
1458: PetscErrorCode ierr;
1461: if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim);
1462: if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits");
1463: /* Find K such that the weights are < 32 digits of precision */
1464: for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) {
1465: wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h)));
1466: }
1467: PetscQuadratureCreate(PETSC_COMM_SELF, q);
1468: PetscQuadratureSetOrder(*q, 2*K+1);
1469: npoints = 2*K-1;
1470: PetscMalloc1(npoints*dim, &x);
1471: PetscMalloc1(npoints, &w);
1472: /* Center term */
1473: x[0] = beta;
1474: w[0] = 0.5*alpha*PETSC_PI;
1475: for (k = 1; k < K; ++k) {
1476: wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1477: xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h));
1478: x[2*k-1] = -alpha*xk+beta;
1479: w[2*k-1] = wk;
1480: x[2*k+0] = alpha*xk+beta;
1481: w[2*k+0] = wk;
1482: }
1483: PetscQuadratureSetData(*q, dim, 1, npoints, x, w);
1484: return(0);
1485: }
1487: PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1488: {
1489: const PetscInt p = 16; /* Digits of precision in the evaluation */
1490: const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */
1491: const PetscReal beta = (b+a)/2.; /* Center of the integration interval */
1492: PetscReal h = 1.0; /* Step size, length between x_k */
1493: PetscInt l = 0; /* Level of refinement, h = 2^{-l} */
1494: PetscReal osum = 0.0; /* Integral on last level */
1495: PetscReal psum = 0.0; /* Integral on the level before the last level */
1496: PetscReal sum; /* Integral on current level */
1497: PetscReal yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */
1498: PetscReal lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */
1499: PetscReal wk; /* Quadrature weight at x_k */
1500: PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */
1501: PetscInt d; /* Digits of precision in the integral */
1504: if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
1505: /* Center term */
1506: func(beta, &lval);
1507: sum = 0.5*alpha*PETSC_PI*lval;
1508: /* */
1509: do {
1510: PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4;
1511: PetscInt k = 1;
1513: ++l;
1514: /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
1515: /* At each level of refinement, h --> h/2 and sum --> sum/2 */
1516: psum = osum;
1517: osum = sum;
1518: h *= 0.5;
1519: sum *= 0.5;
1520: do {
1521: wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1522: yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1523: lx = -alpha*(1.0 - yk)+beta;
1524: rx = alpha*(1.0 - yk)+beta;
1525: func(lx, &lval);
1526: func(rx, &rval);
1527: lterm = alpha*wk*lval;
1528: maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm);
1529: sum += lterm;
1530: rterm = alpha*wk*rval;
1531: maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm);
1532: sum += rterm;
1533: ++k;
1534: /* Only need to evaluate every other point on refined levels */
1535: if (l != 1) ++k;
1536: } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */
1538: d1 = PetscLog10Real(PetscAbsReal(sum - osum));
1539: d2 = PetscLog10Real(PetscAbsReal(sum - psum));
1540: d3 = PetscLog10Real(maxTerm) - p;
1541: if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0;
1542: else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)));
1543: d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
1544: } while (d < digits && l < 12);
1545: *sol = sum;
1547: return(0);
1548: }
1550: #if defined(PETSC_HAVE_MPFR)
1551: PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1552: {
1553: const PetscInt safetyFactor = 2; /* Calculate abcissa until 2*p digits */
1554: PetscInt l = 0; /* Level of refinement, h = 2^{-l} */
1555: mpfr_t alpha; /* Half-width of the integration interval */
1556: mpfr_t beta; /* Center of the integration interval */
1557: mpfr_t h; /* Step size, length between x_k */
1558: mpfr_t osum; /* Integral on last level */
1559: mpfr_t psum; /* Integral on the level before the last level */
1560: mpfr_t sum; /* Integral on current level */
1561: mpfr_t yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */
1562: mpfr_t lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */
1563: mpfr_t wk; /* Quadrature weight at x_k */
1564: PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */
1565: PetscInt d; /* Digits of precision in the integral */
1566: mpfr_t pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp;
1569: if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
1570: /* Create high precision storage */
1571: 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);
1572: /* Initialization */
1573: mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN);
1574: mpfr_set_d(beta, 0.5*(b+a), MPFR_RNDN);
1575: mpfr_set_d(osum, 0.0, MPFR_RNDN);
1576: mpfr_set_d(psum, 0.0, MPFR_RNDN);
1577: mpfr_set_d(h, 1.0, MPFR_RNDN);
1578: mpfr_const_pi(pi2, MPFR_RNDN);
1579: mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN);
1580: /* Center term */
1581: func(0.5*(b+a), &lval);
1582: mpfr_set(sum, pi2, MPFR_RNDN);
1583: mpfr_mul(sum, sum, alpha, MPFR_RNDN);
1584: mpfr_mul_d(sum, sum, lval, MPFR_RNDN);
1585: /* */
1586: do {
1587: PetscReal d1, d2, d3, d4;
1588: PetscInt k = 1;
1590: ++l;
1591: mpfr_set_d(maxTerm, 0.0, MPFR_RNDN);
1592: /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
1593: /* At each level of refinement, h --> h/2 and sum --> sum/2 */
1594: mpfr_set(psum, osum, MPFR_RNDN);
1595: mpfr_set(osum, sum, MPFR_RNDN);
1596: mpfr_mul_d(h, h, 0.5, MPFR_RNDN);
1597: mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN);
1598: do {
1599: mpfr_set_si(kh, k, MPFR_RNDN);
1600: mpfr_mul(kh, kh, h, MPFR_RNDN);
1601: /* Weight */
1602: mpfr_set(wk, h, MPFR_RNDN);
1603: mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN);
1604: mpfr_mul(msinh, msinh, pi2, MPFR_RNDN);
1605: mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN);
1606: mpfr_cosh(tmp, msinh, MPFR_RNDN);
1607: mpfr_sqr(tmp, tmp, MPFR_RNDN);
1608: mpfr_mul(wk, wk, mcosh, MPFR_RNDN);
1609: mpfr_div(wk, wk, tmp, MPFR_RNDN);
1610: /* Abscissa */
1611: mpfr_set_d(yk, 1.0, MPFR_RNDZ);
1612: mpfr_cosh(tmp, msinh, MPFR_RNDN);
1613: mpfr_div(yk, yk, tmp, MPFR_RNDZ);
1614: mpfr_exp(tmp, msinh, MPFR_RNDN);
1615: mpfr_div(yk, yk, tmp, MPFR_RNDZ);
1616: /* Quadrature points */
1617: mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ);
1618: mpfr_mul(lx, lx, alpha, MPFR_RNDU);
1619: mpfr_add(lx, lx, beta, MPFR_RNDU);
1620: mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ);
1621: mpfr_mul(rx, rx, alpha, MPFR_RNDD);
1622: mpfr_add(rx, rx, beta, MPFR_RNDD);
1623: /* Evaluation */
1624: func(mpfr_get_d(lx, MPFR_RNDU), &lval);
1625: func(mpfr_get_d(rx, MPFR_RNDD), &rval);
1626: /* Update */
1627: mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
1628: mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN);
1629: mpfr_add(sum, sum, tmp, MPFR_RNDN);
1630: mpfr_abs(tmp, tmp, MPFR_RNDN);
1631: mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
1632: mpfr_set(curTerm, tmp, MPFR_RNDN);
1633: mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
1634: mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN);
1635: mpfr_add(sum, sum, tmp, MPFR_RNDN);
1636: mpfr_abs(tmp, tmp, MPFR_RNDN);
1637: mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
1638: mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN);
1639: ++k;
1640: /* Only need to evaluate every other point on refined levels */
1641: if (l != 1) ++k;
1642: mpfr_log10(tmp, wk, MPFR_RNDN);
1643: mpfr_abs(tmp, tmp, MPFR_RNDN);
1644: } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */
1645: mpfr_sub(tmp, sum, osum, MPFR_RNDN);
1646: mpfr_abs(tmp, tmp, MPFR_RNDN);
1647: mpfr_log10(tmp, tmp, MPFR_RNDN);
1648: d1 = mpfr_get_d(tmp, MPFR_RNDN);
1649: mpfr_sub(tmp, sum, psum, MPFR_RNDN);
1650: mpfr_abs(tmp, tmp, MPFR_RNDN);
1651: mpfr_log10(tmp, tmp, MPFR_RNDN);
1652: d2 = mpfr_get_d(tmp, MPFR_RNDN);
1653: mpfr_log10(tmp, maxTerm, MPFR_RNDN);
1654: d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits;
1655: mpfr_log10(tmp, curTerm, MPFR_RNDN);
1656: d4 = mpfr_get_d(tmp, MPFR_RNDN);
1657: d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
1658: } while (d < digits && l < 8);
1659: *sol = mpfr_get_d(sum, MPFR_RNDN);
1660: /* Cleanup */
1661: mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
1662: return(0);
1663: }
1664: #else
1666: PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1667: {
1668: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp");
1669: }
1670: #endif
1672: /* Overwrites A. Can only handle full-rank problems with m>=n
1673: * A in column-major format
1674: * Ainv in row-major format
1675: * tau has length m
1676: * worksize must be >= max(1,n)
1677: */
1678: static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
1679: {
1681: PetscBLASInt M,N,K,lda,ldb,ldwork,info;
1682: PetscScalar *A,*Ainv,*R,*Q,Alpha;
1685: #if defined(PETSC_USE_COMPLEX)
1686: {
1687: PetscInt i,j;
1688: PetscMalloc2(m*n,&A,m*n,&Ainv);
1689: for (j=0; j<n; j++) {
1690: for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
1691: }
1692: mstride = m;
1693: }
1694: #else
1695: A = A_in;
1696: Ainv = Ainv_out;
1697: #endif
1699: PetscBLASIntCast(m,&M);
1700: PetscBLASIntCast(n,&N);
1701: PetscBLASIntCast(mstride,&lda);
1702: PetscBLASIntCast(worksize,&ldwork);
1703: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1704: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
1705: PetscFPTrapPop();
1706: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
1707: R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
1709: /* Extract an explicit representation of Q */
1710: Q = Ainv;
1711: PetscArraycpy(Q,A,mstride*n);
1712: K = N; /* full rank */
1713: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
1714: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
1716: /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
1717: Alpha = 1.0;
1718: ldb = lda;
1719: PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
1720: /* Ainv is Q, overwritten with inverse */
1722: #if defined(PETSC_USE_COMPLEX)
1723: {
1724: PetscInt i;
1725: for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
1726: PetscFree2(A,Ainv);
1727: }
1728: #endif
1729: return(0);
1730: }
1732: /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
1733: static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
1734: {
1736: PetscReal *Bv;
1737: PetscInt i,j;
1740: PetscMalloc1((ninterval+1)*ndegree,&Bv);
1741: /* Point evaluation of L_p on all the source vertices */
1742: PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);
1743: /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
1744: for (i=0; i<ninterval; i++) {
1745: for (j=0; j<ndegree; j++) {
1746: if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1747: else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1748: }
1749: }
1750: PetscFree(Bv);
1751: return(0);
1752: }
1754: /*@
1755: PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
1757: Not Collective
1759: Input Arguments:
1760: + degree - degree of reconstruction polynomial
1761: . nsource - number of source intervals
1762: . sourcex - sorted coordinates of source cell boundaries (length nsource+1)
1763: . ntarget - number of target intervals
1764: - targetx - sorted coordinates of target cell boundaries (length ntarget+1)
1766: Output Arguments:
1767: . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
1769: Level: advanced
1771: .seealso: PetscDTLegendreEval()
1772: @*/
1773: PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
1774: {
1776: PetscInt i,j,k,*bdegrees,worksize;
1777: PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
1778: PetscScalar *tau,*work;
1784: 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);
1785: #if defined(PETSC_USE_DEBUG)
1786: for (i=0; i<nsource; i++) {
1787: 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]);
1788: }
1789: for (i=0; i<ntarget; i++) {
1790: 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]);
1791: }
1792: #endif
1793: xmin = PetscMin(sourcex[0],targetx[0]);
1794: xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
1795: center = (xmin + xmax)/2;
1796: hscale = (xmax - xmin)/2;
1797: worksize = nsource;
1798: PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);
1799: PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);
1800: for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
1801: for (i=0; i<=degree; i++) bdegrees[i] = i+1;
1802: PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);
1803: PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);
1804: for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
1805: PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);
1806: for (i=0; i<ntarget; i++) {
1807: PetscReal rowsum = 0;
1808: for (j=0; j<nsource; j++) {
1809: PetscReal sum = 0;
1810: for (k=0; k<degree+1; k++) {
1811: sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
1812: }
1813: R[i*nsource+j] = sum;
1814: rowsum += sum;
1815: }
1816: for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
1817: }
1818: PetscFree4(bdegrees,sourcey,Bsource,work);
1819: PetscFree4(tau,Bsinv,targety,Btarget);
1820: return(0);
1821: }
1823: /*@C
1824: PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points
1826: Not Collective
1828: Input Parameter:
1829: + n - the number of GLL nodes
1830: . nodes - the GLL nodes
1831: . weights - the GLL weights
1832: - f - the function values at the nodes
1834: Output Parameter:
1835: . in - the value of the integral
1837: Level: beginner
1839: .seealso: PetscDTGaussLobattoLegendreQuadrature()
1841: @*/
1842: PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n,PetscReal *nodes,PetscReal *weights,const PetscReal *f,PetscReal *in)
1843: {
1844: PetscInt i;
1847: *in = 0.;
1848: for (i=0; i<n; i++) {
1849: *in += f[i]*f[i]*weights[i];
1850: }
1851: return(0);
1852: }
1854: /*@C
1855: PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element
1857: Not Collective
1859: Input Parameter:
1860: + n - the number of GLL nodes
1861: . nodes - the GLL nodes
1862: - weights - the GLL weights
1864: Output Parameter:
1865: . A - the stiffness element
1867: Level: beginner
1869: Notes:
1870: Destroy this with PetscGaussLobattoLegendreElementLaplacianDestroy()
1872: You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented (the array is symmetric)
1874: .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy()
1876: @*/
1877: PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1878: {
1879: PetscReal **A;
1880: PetscErrorCode ierr;
1881: const PetscReal *gllnodes = nodes;
1882: const PetscInt p = n-1;
1883: PetscReal z0,z1,z2 = -1,x,Lpj,Lpr;
1884: PetscInt i,j,nn,r;
1887: PetscMalloc1(n,&A);
1888: PetscMalloc1(n*n,&A[0]);
1889: for (i=1; i<n; i++) A[i] = A[i-1]+n;
1891: for (j=1; j<p; j++) {
1892: x = gllnodes[j];
1893: z0 = 1.;
1894: z1 = x;
1895: for (nn=1; nn<p; nn++) {
1896: z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1897: z0 = z1;
1898: z1 = z2;
1899: }
1900: Lpj=z2;
1901: for (r=1; r<p; r++) {
1902: if (r == j) {
1903: A[j][j]=2./(3.*(1.-gllnodes[j]*gllnodes[j])*Lpj*Lpj);
1904: } else {
1905: x = gllnodes[r];
1906: z0 = 1.;
1907: z1 = x;
1908: for (nn=1; nn<p; nn++) {
1909: z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1910: z0 = z1;
1911: z1 = z2;
1912: }
1913: Lpr = z2;
1914: A[r][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*Lpr*(gllnodes[j]-gllnodes[r])*(gllnodes[j]-gllnodes[r]));
1915: }
1916: }
1917: }
1918: for (j=1; j<p+1; j++) {
1919: x = gllnodes[j];
1920: z0 = 1.;
1921: z1 = x;
1922: for (nn=1; nn<p; nn++) {
1923: z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1924: z0 = z1;
1925: z1 = z2;
1926: }
1927: Lpj = z2;
1928: A[j][0] = 4.*PetscPowRealInt(-1.,p)/(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.+gllnodes[j])*(1.+gllnodes[j]));
1929: A[0][j] = A[j][0];
1930: }
1931: for (j=0; j<p; j++) {
1932: x = gllnodes[j];
1933: z0 = 1.;
1934: z1 = x;
1935: for (nn=1; nn<p; nn++) {
1936: z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1937: z0 = z1;
1938: z1 = z2;
1939: }
1940: Lpj=z2;
1942: A[p][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.-gllnodes[j])*(1.-gllnodes[j]));
1943: A[j][p] = A[p][j];
1944: }
1945: A[0][0]=0.5+(((PetscReal)p)*(((PetscReal)p)+1.)-2.)/6.;
1946: A[p][p]=A[0][0];
1947: *AA = A;
1948: return(0);
1949: }
1951: /*@C
1952: PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element
1954: Not Collective
1956: Input Parameter:
1957: + n - the number of GLL nodes
1958: . nodes - the GLL nodes
1959: . weights - the GLL weightss
1960: - A - the stiffness element
1962: Level: beginner
1964: .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate()
1966: @*/
1967: PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1968: {
1972: PetscFree((*AA)[0]);
1973: PetscFree(*AA);
1974: *AA = NULL;
1975: return(0);
1976: }
1978: /*@C
1979: PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element
1981: Not Collective
1983: Input Parameter:
1984: + n - the number of GLL nodes
1985: . nodes - the GLL nodes
1986: . weights - the GLL weights
1988: Output Parameter:
1989: . AA - the stiffness element
1990: - AAT - the transpose of AA (pass in NULL if you do not need this array)
1992: Level: beginner
1994: Notes:
1995: Destroy this with PetscGaussLobattoLegendreElementGradientDestroy()
1997: You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented
1999: .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy()
2001: @*/
2002: PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT)
2003: {
2004: PetscReal **A, **AT = NULL;
2005: PetscErrorCode ierr;
2006: const PetscReal *gllnodes = nodes;
2007: const PetscInt p = n-1;
2008: PetscReal Li, Lj,d0;
2009: PetscInt i,j;
2012: PetscMalloc1(n,&A);
2013: PetscMalloc1(n*n,&A[0]);
2014: for (i=1; i<n; i++) A[i] = A[i-1]+n;
2016: if (AAT) {
2017: PetscMalloc1(n,&AT);
2018: PetscMalloc1(n*n,&AT[0]);
2019: for (i=1; i<n; i++) AT[i] = AT[i-1]+n;
2020: }
2022: if (n==1) {A[0][0] = 0.;}
2023: d0 = (PetscReal)p*((PetscReal)p+1.)/4.;
2024: for (i=0; i<n; i++) {
2025: for (j=0; j<n; j++) {
2026: A[i][j] = 0.;
2027: PetscDTComputeJacobi(0., 0., p, gllnodes[i], &Li);
2028: PetscDTComputeJacobi(0., 0., p, gllnodes[j], &Lj);
2029: if (i!=j) A[i][j] = Li/(Lj*(gllnodes[i]-gllnodes[j]));
2030: if ((j==i) && (i==0)) A[i][j] = -d0;
2031: if (j==i && i==p) A[i][j] = d0;
2032: if (AT) AT[j][i] = A[i][j];
2033: }
2034: }
2035: if (AAT) *AAT = AT;
2036: *AA = A;
2037: return(0);
2038: }
2040: /*@C
2041: PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with PetscGaussLobattoLegendreElementGradientCreate()
2043: Not Collective
2045: Input Parameter:
2046: + n - the number of GLL nodes
2047: . nodes - the GLL nodes
2048: . weights - the GLL weights
2049: . AA - the stiffness element
2050: - AAT - the transpose of the element
2052: Level: beginner
2054: .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionCreate()
2056: @*/
2057: PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT)
2058: {
2062: PetscFree((*AA)[0]);
2063: PetscFree(*AA);
2064: *AA = NULL;
2065: if (*AAT) {
2066: PetscFree((*AAT)[0]);
2067: PetscFree(*AAT);
2068: *AAT = NULL;
2069: }
2070: return(0);
2071: }
2073: /*@C
2074: PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element
2076: Not Collective
2078: Input Parameter:
2079: + n - the number of GLL nodes
2080: . nodes - the GLL nodes
2081: - weights - the GLL weightss
2083: Output Parameter:
2084: . AA - the stiffness element
2086: Level: beginner
2088: Notes:
2089: Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy()
2091: This is the same as the Gradient operator multiplied by the diagonal mass matrix
2093: You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented
2095: .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy()
2097: @*/
2098: PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
2099: {
2100: PetscReal **D;
2101: PetscErrorCode ierr;
2102: const PetscReal *gllweights = weights;
2103: const PetscInt glln = n;
2104: PetscInt i,j;
2107: PetscGaussLobattoLegendreElementGradientCreate(n,nodes,weights,&D,NULL);
2108: for (i=0; i<glln; i++){
2109: for (j=0; j<glln; j++) {
2110: D[i][j] = gllweights[i]*D[i][j];
2111: }
2112: }
2113: *AA = D;
2114: return(0);
2115: }
2117: /*@C
2118: PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element
2120: Not Collective
2122: Input Parameter:
2123: + n - the number of GLL nodes
2124: . nodes - the GLL nodes
2125: . weights - the GLL weights
2126: - A - advection
2128: Level: beginner
2130: .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementAdvectionCreate()
2132: @*/
2133: PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
2134: {
2138: PetscFree((*AA)[0]);
2139: PetscFree(*AA);
2140: *AA = NULL;
2141: return(0);
2142: }
2144: PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
2145: {
2146: PetscReal **A;
2147: PetscErrorCode ierr;
2148: const PetscReal *gllweights = weights;
2149: const PetscInt glln = n;
2150: PetscInt i,j;
2153: PetscMalloc1(glln,&A);
2154: PetscMalloc1(glln*glln,&A[0]);
2155: for (i=1; i<glln; i++) A[i] = A[i-1]+glln;
2156: if (glln==1) {A[0][0] = 0.;}
2157: for (i=0; i<glln; i++) {
2158: for (j=0; j<glln; j++) {
2159: A[i][j] = 0.;
2160: if (j==i) A[i][j] = gllweights[i];
2161: }
2162: }
2163: *AA = A;
2164: return(0);
2165: }
2167: PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
2168: {
2172: PetscFree((*AA)[0]);
2173: PetscFree(*AA);
2174: *AA = NULL;
2175: return(0);
2176: }
2178: /*@
2179: PetscDTIndexToBary - convert an index into a barycentric coordinate.
2181: Input Parameters:
2182: + len - the desired length of the barycentric tuple (usually 1 more than the dimension it represents, so a barycentric coordinate in a triangle has length 3)
2183: . sum - the value that the sum of the barycentric coordinates (which will be non-negative integers) should sum to
2184: - index - the index to convert: should be >= 0 and < Binomial(len - 1 + sum, sum)
2186: Output Parameter:
2187: . coord - will be filled with the barycentric coordinate
2189: Level: beginner
2191: Note: the indices map to barycentric coordinates in lexicographic order, where the first index is the
2192: least significant and the last index is the most significant.
2194: .seealso: PetscDTBaryToIndex
2195: @*/
2196: PetscErrorCode PetscDTIndexToBary(PetscInt len, PetscInt sum, PetscInt index, PetscInt coord[])
2197: {
2198: PetscInt c, d, s, total, subtotal, nexttotal;
2201: if (len < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative");
2202: if (index < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index must be non-negative");
2203: if (!len) {
2204: if (!sum && !index) return(0);
2205: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid index or sum for length 0 barycentric coordinate");
2206: }
2207: for (c = 1, total = 1; c <= len; c++) {
2208: /* total is the number of ways to have a tuple of length c with sum */
2209: if (index < total) break;
2210: total = (total * (sum + c)) / c;
2211: }
2212: if (c > len) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index out of range");
2213: for (d = c; d < len; d++) coord[d] = 0;
2214: for (s = 0, subtotal = 1, nexttotal = 1; c > 0;) {
2215: /* subtotal is the number of ways to have a tuple of length c with sum s */
2216: /* nexttotal is the number of ways to have a tuple of length c-1 with sum s */
2217: if ((index + subtotal) >= total) {
2218: coord[--c] = sum - s;
2219: index -= (total - subtotal);
2220: sum = s;
2221: total = nexttotal;
2222: subtotal = 1;
2223: nexttotal = 1;
2224: s = 0;
2225: } else {
2226: subtotal = (subtotal * (c + s)) / (s + 1);
2227: nexttotal = (nexttotal * (c - 1 + s)) / (s + 1);
2228: s++;
2229: }
2230: }
2231: return(0);
2232: }
2234: /*@
2235: PetscDTBaryToIndex - convert a barycentric coordinate to an index
2237: Input Parameters:
2238: + len - the desired length of the barycentric tuple (usually 1 more than the dimension it represents, so a barycentric coordinate in a triangle has length 3)
2239: . sum - the value that the sum of the barycentric coordinates (which will be non-negative integers) should sum to
2240: - coord - a barycentric coordinate with the given length and sum
2242: Output Parameter:
2243: . index - the unique index for the coordinate, >= 0 and < Binomial(len - 1 + sum, sum)
2245: Level: beginner
2247: Note: the indices map to barycentric coordinates in lexicographic order, where the first index is the
2248: least significant and the last index is the most significant.
2250: .seealso: PetscDTIndexToBary
2251: @*/
2252: PetscErrorCode PetscDTBaryToIndex(PetscInt len, PetscInt sum, const PetscInt coord[], PetscInt *index)
2253: {
2254: PetscInt c;
2255: PetscInt i;
2256: PetscInt total;
2259: if (len < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative");
2260: if (!len) {
2261: if (!sum) {
2262: *index = 0;
2263: return(0);
2264: }
2265: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid index or sum for length 0 barycentric coordinate");
2266: }
2267: for (c = 1, total = 1; c < len; c++) total = (total * (sum + c)) / c;
2268: i = total - 1;
2269: c = len - 1;
2270: sum -= coord[c];
2271: while (sum > 0) {
2272: PetscInt subtotal;
2273: PetscInt s;
2275: for (s = 1, subtotal = 1; s < sum; s++) subtotal = (subtotal * (c + s)) / s;
2276: i -= subtotal;
2277: sum -= coord[--c];
2278: }
2279: *index = i;
2280: return(0);
2281: }