Actual source code: ex1.c
petsc-3.14.0 2020-09-29
1: static char help[] = "Tests 1D discretization tools.\n\n";
3: #include <petscdt.h>
4: #include <petscviewer.h>
5: #include <petsc/private/petscimpl.h>
6: #include <petsc/private/dtimpl.h>
8: static PetscErrorCode CheckPoints(const char *name,PetscInt npoints,const PetscReal *points,PetscInt ndegrees,const PetscInt *degrees)
9: {
11: PetscReal *B,*D,*D2;
12: PetscInt i,j;
15: PetscMalloc3(npoints*ndegrees,&B,npoints*ndegrees,&D,npoints*ndegrees,&D2);
16: PetscDTLegendreEval(npoints,points,ndegrees,degrees,B,D,D2);
17: PetscPrintf(PETSC_COMM_WORLD,"%s\n",name);
18: for (i=0; i<npoints; i++) {
19: for (j=0; j<ndegrees; j++) {
20: PetscReal b,d,d2;
21: b = B[i*ndegrees+j];
22: d = D[i*ndegrees+j];
23: d2 = D2[i*ndegrees+j];
24: if (PetscAbsReal(b) < PETSC_SMALL) b = 0;
25: if (PetscAbsReal(d) < PETSC_SMALL) d = 0;
26: if (PetscAbsReal(d2) < PETSC_SMALL) d2 = 0;
27: PetscPrintf(PETSC_COMM_WORLD,"degree %D at %12.4g: B=%12.4g D=%12.4g D2=%12.4g\n",degrees[j],(double)points[i],(double)b,(double)d,(double)d2);
28: }
29: }
30: PetscFree3(B,D,D2);
31: return(0);
32: }
34: typedef PetscErrorCode(*quadratureFunc)(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal[],PetscReal[]);
36: static PetscErrorCode CheckQuadrature_Basics(PetscInt npoints, PetscReal alpha, PetscReal beta, const PetscReal x[], const PetscReal w[])
37: {
38: PetscInt i;
41: for (i = 1; i < npoints; i++) {
42: if (x[i] <= x[i-1]) SETERRQ6(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Quadrature points not monotonically increasing, %D points, alpha = %g, beta = %g, i = %D, x[i] = %g, x[i-1] = %g\n",npoints, (double) alpha, (double) beta, i, x[i], x[i-1]);
43: }
44: for (i = 0; i < npoints; i++) {
45: if (w[i] <= 0.) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Quadrature weight not positive, %D points, alpha = %g, beta = %g, i = %D, w[i] = %g\n",npoints, (double) alpha, (double) beta, i, w[i]);
46: }
47: return(0);
48: }
50: static PetscErrorCode CheckQuadrature(PetscInt npoints, PetscReal alpha, PetscReal beta, const PetscReal x[], const PetscReal w[], PetscInt nexact)
51: {
52: PetscInt i, j, k;
53: PetscReal *Pi, *Pj;
54: PetscReal eps;
58: eps = PETSC_SMALL;
59: PetscMalloc2(npoints, &Pi, npoints, &Pj);
60: for (i = 0; i <= nexact; i++) {
61: PetscDTJacobiEval(npoints, alpha, beta, x, 1, &i, Pi, NULL, NULL);
62: for (j = i; j <= nexact - i; j++) {
63: PetscReal I_quad = 0.;
64: PetscReal I_exact = 0.;
65: PetscReal err, tol;
66: PetscDTJacobiEval(npoints, alpha, beta, x, 1, &j, Pj, NULL, NULL);
68: tol = eps;
69: if (i == j) {
70: PetscReal norm, norm2diff;
72: I_exact = PetscPowReal(2.0, alpha + beta + 1.) / (2.*i + alpha + beta + 1.);
73: #if defined(PETSC_HAVE_LGAMMA)
74: I_exact *= PetscExpReal(PetscLGamma(i + alpha + 1.) + PetscLGamma(i + beta + 1.) - (PetscLGamma(i + alpha + beta + 1.) + PetscLGamma(i + 1.)));
75: #else
76: {
77: PetscInt ibeta = (PetscInt) beta;
79: if ((PetscReal) ibeta != beta) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable.");
80: for (k = 0; k < ibeta; k++) I_exact *= (i + 1. + k) / (i + alpha + 1. + k);
81: }
82: #endif
84: PetscDTJacobiNorm(alpha, beta, i, &norm);
85: norm2diff = PetscAbsReal(norm*norm - I_exact);
86: if (norm2diff > eps * I_exact) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Jacobi norm error %g\n", (double) norm2diff);
88: tol = eps * I_exact;
89: }
90: for (k = 0; k < npoints; k++) I_quad += w[k] * (Pi[k] * Pj[k]);
91: err = PetscAbsReal(I_exact - I_quad);
92: PetscInfo7(NULL,"npoints %D, alpha %g, beta %g, i %D, j %D, exact %g, err %g\n", npoints, (double) alpha, (double) beta, i, j, (double) I_exact, (double) err);
93: if (err > tol) SETERRQ7(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrectly integrated P_%D * P_%D using %D point rule with alpha = %g, beta = %g: exact %g, err %g\n", i, j, npoints, (double) alpha, (double) beta, (double) I_exact, (double) err);
94: }
95: }
96: PetscFree2(Pi, Pj);
97: return(0);
98: }
100: static PetscErrorCode CheckJacobiQuadrature(PetscInt npoints, PetscReal alpha, PetscReal beta, quadratureFunc func, PetscInt nexact)
101: {
102: PetscReal *x, *w;
106: PetscMalloc2(npoints, &x, npoints, &w);
107: (*func)(npoints, -1., 1., alpha, beta, x, w);
108: CheckQuadrature_Basics(npoints, alpha, beta, x, w);
109: CheckQuadrature(npoints, alpha, beta, x, w, nexact);
110: #if defined(PETSCDTGAUSSIANQUADRATURE_EIG)
111: /* compare methods of computing quadrature */
112: PetscDTGaussQuadratureNewton_Internal = (PetscBool) !PetscDTGaussQuadratureNewton_Internal;
113: {
114: PetscReal *x2, *w2;
115: PetscReal eps;
116: PetscInt i;
118: eps = PETSC_SMALL;
119: PetscMalloc2(npoints, &x2, npoints, &w2);
120: (*func)(npoints, -1., 1., alpha, beta, x2, w2);
121: CheckQuadrature_Basics(npoints, alpha, beta, x2, w2);
122: CheckQuadrature(npoints, alpha, beta, x2, w2, nexact);
123: for (i = 0; i < npoints; i++) {
124: PetscReal xdiff, xtol, wdiff, wtol;
126: xdiff = PetscAbsReal(x[i] - x2[i]);
127: wdiff = PetscAbsReal(w[i] - w2[i]);
128: xtol = eps * (1. + PetscMin(PetscAbsReal(x[i]),1. - PetscAbsReal(x[i])));
129: wtol = eps * (1. + w[i]);
130: PetscInfo6(NULL,"npoints %D, alpha %g, beta %g, i %D, xdiff/xtol %g, wdiff/wtol %g\n", npoints, (double) alpha, (double) beta, i, (double) xdiff/xtol, (double) wdiff/wtol);
131: if (xdiff > xtol) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Mismatch quadrature point: %D points, alpha = %g, beta = %g, i = %D, xdiff = %g\n", npoints, (double) alpha, (double) beta, i, (double) xdiff);
132: if (wdiff > wtol) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Mismatch quadrature weight: %D points, alpha = %g, beta = %g, i = %D, wdiff = %g\n", npoints, (double) alpha, (double) beta, i, (double) wdiff);
133: }
134: PetscFree2(x2, w2);
135: }
136: /* restore */
137: PetscDTGaussQuadratureNewton_Internal = (PetscBool) !PetscDTGaussQuadratureNewton_Internal;
138: #endif
139: PetscFree2(x, w);
140: return(0);
141: }
143: int main(int argc,char **argv)
144: {
146: PetscInt degrees[1000],ndegrees,npoints,two;
147: PetscReal points[1000],weights[1000],interval[2];
148: PetscInt minpoints, maxpoints;
149: PetscBool flg;
151: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
152: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Discretization tools test options",NULL);
153: {
154: ndegrees = 1000;
155: degrees[0] = 0;
156: degrees[1] = 1;
157: degrees[2] = 2;
158: PetscOptionsIntArray("-degrees","list of degrees to evaluate","",degrees,&ndegrees,&flg);
160: if (!flg) ndegrees = 3;
161: npoints = 1000;
162: points[0] = 0.0;
163: points[1] = -0.5;
164: points[2] = 1.0;
165: PetscOptionsRealArray("-points","list of points at which to evaluate","",points,&npoints,&flg);
167: if (!flg) npoints = 3;
168: two = 2;
169: interval[0] = -1.;
170: interval[1] = 1.;
171: PetscOptionsRealArray("-interval","interval on which to construct quadrature","",interval,&two,NULL);
173: minpoints = 1;
174: PetscOptionsInt("-minpoints","minimum points for thorough Gauss-Jacobi quadrature tests","",minpoints,&minpoints,NULL);
175: maxpoints = 30;
176: #if defined(PETSC_USE_REAL_SINGLE)
177: maxpoints = 5;
178: #elif defined(PETSC_USE_REAL___FLOAT128)
179: maxpoints = 20; /* just to make test faster */
180: #endif
181: PetscOptionsInt("-maxpoints","maximum points for thorough Gauss-Jacobi quadrature tests","",maxpoints,&maxpoints,NULL);
182: }
183: PetscOptionsEnd();
184: CheckPoints("User-provided points",npoints,points,ndegrees,degrees);
186: PetscDTGaussQuadrature(npoints,interval[0],interval[1],points,weights);
187: PetscPrintf(PETSC_COMM_WORLD,"Quadrature weights\n");
188: PetscRealView(npoints,weights,PETSC_VIEWER_STDOUT_WORLD);
189: {
190: PetscReal a = interval[0],b = interval[1],zeroth,first,second;
191: PetscInt i;
192: zeroth = b - a;
193: first = (b*b - a*a)/2;
194: second = (b*b*b - a*a*a)/3;
195: for (i=0; i<npoints; i++) {
196: zeroth -= weights[i];
197: first -= weights[i] * points[i];
198: second -= weights[i] * PetscSqr(points[i]);
199: }
200: if (PetscAbs(zeroth) < 1e-10) zeroth = 0.;
201: if (PetscAbs(first) < 1e-10) first = 0.;
202: if (PetscAbs(second) < 1e-10) second = 0.;
203: PetscPrintf(PETSC_COMM_WORLD,"Moment error: zeroth=%g, first=%g, second=%g\n",(double)(-zeroth),(double)(-first),(double)(-second));
204: }
205: CheckPoints("Gauss points",npoints,points,ndegrees,degrees);
206: {
207: PetscInt i;
209: for (i = minpoints; i <= maxpoints; i++) {
210: PetscReal a1, b1, a2, b2;
212: #if defined(PETSC_HAVE_LGAMMA)
213: a1 = -0.6;
214: b1 = 1.1;
215: a2 = 2.2;
216: b2 = -0.6;
217: #else
218: a1 = 0.;
219: b1 = 1.;
220: a2 = 2.;
221: b2 = 0.;
222: #endif
223: CheckJacobiQuadrature(i, 0., 0., PetscDTGaussJacobiQuadrature, 2*i-1);
224: CheckJacobiQuadrature(i, a1, b1, PetscDTGaussJacobiQuadrature, 2*i-1);
225: CheckJacobiQuadrature(i, a2, b2, PetscDTGaussJacobiQuadrature, 2*i-1);
226: if (i >= 2) {
227: CheckJacobiQuadrature(i, 0., 0., PetscDTGaussLobattoJacobiQuadrature, 2*i-3);
228: CheckJacobiQuadrature(i, a1, b1, PetscDTGaussLobattoJacobiQuadrature, 2*i-3);
229: CheckJacobiQuadrature(i, a2, b2, PetscDTGaussLobattoJacobiQuadrature, 2*i-3);
230: }
231: }
232: }
233: PetscFinalize();
234: return ierr;
235: }
237: /*TEST
238: test:
239: suffix: 1
240: args: -degrees 1,2,3,4,5 -points 0,.2,-.5,.8,.9,1 -interval -.5,1
241: TEST*/