static char help[] = "Tests 1D discretization tools.\n\n"; #include #include #include #include static PetscErrorCode CheckPoints(const char *name,PetscInt npoints,const PetscReal *points,PetscInt ndegrees,const PetscInt *degrees) { PetscReal *B,*D,*D2; PetscInt i,j; PetscFunctionBegin; PetscCall(PetscMalloc3(npoints*ndegrees,&B,npoints*ndegrees,&D,npoints*ndegrees,&D2)); PetscCall(PetscDTLegendreEval(npoints,points,ndegrees,degrees,B,D,D2)); PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s\n",name)); for (i=0; i eps * I_exact,PETSC_COMM_SELF,PETSC_ERR_PLIB, "Jacobi norm error %g", (double) norm2diff); tol = eps * I_exact; } for (k = 0; k < npoints; k++) I_quad += w[k] * (Pi[k] * Pj[k]); err = PetscAbsReal(I_exact - I_quad); PetscCall(PetscInfo(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)); PetscCheckFalse(err > tol,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", i, j, npoints, (double) alpha, (double) beta, (double) I_exact, (double) err); } } PetscCall(PetscFree2(Pi, Pj)); PetscFunctionReturn(0); } static PetscErrorCode CheckJacobiQuadrature(PetscInt npoints, PetscReal alpha, PetscReal beta, quadratureFunc func, PetscInt nexact) { PetscReal *x, *w; PetscFunctionBegin; PetscCall(PetscMalloc2(npoints, &x, npoints, &w)); PetscCall((*func)(npoints, -1., 1., alpha, beta, x, w)); PetscCall(CheckQuadrature_Basics(npoints, alpha, beta, x, w)); PetscCall(CheckQuadrature(npoints, alpha, beta, x, w, nexact)); #if defined(PETSCDTGAUSSIANQUADRATURE_EIG) /* compare methods of computing quadrature */ PetscDTGaussQuadratureNewton_Internal = (PetscBool) !PetscDTGaussQuadratureNewton_Internal; { PetscReal *x2, *w2; PetscReal eps; PetscInt i; eps = PETSC_SMALL; PetscCall(PetscMalloc2(npoints, &x2, npoints, &w2)); PetscCall((*func)(npoints, -1., 1., alpha, beta, x2, w2)); PetscCall(CheckQuadrature_Basics(npoints, alpha, beta, x2, w2)); PetscCall(CheckQuadrature(npoints, alpha, beta, x2, w2, nexact)); for (i = 0; i < npoints; i++) { PetscReal xdiff, xtol, wdiff, wtol; xdiff = PetscAbsReal(x[i] - x2[i]); wdiff = PetscAbsReal(w[i] - w2[i]); xtol = eps * (1. + PetscMin(PetscAbsReal(x[i]),1. - PetscAbsReal(x[i]))); wtol = eps * (1. + w[i]); PetscCall(PetscInfo(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)); PetscCheckFalse(xdiff > xtol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Mismatch quadrature point: %D points, alpha = %g, beta = %g, i = %D, xdiff = %g", npoints, (double) alpha, (double) beta, i, (double) xdiff); PetscCheckFalse(wdiff > wtol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Mismatch quadrature weight: %D points, alpha = %g, beta = %g, i = %D, wdiff = %g", npoints, (double) alpha, (double) beta, i, (double) wdiff); } PetscCall(PetscFree2(x2, w2)); } /* restore */ PetscDTGaussQuadratureNewton_Internal = (PetscBool) !PetscDTGaussQuadratureNewton_Internal; #endif PetscCall(PetscFree2(x, w)); PetscFunctionReturn(0); } int main(int argc,char **argv) { PetscErrorCode ierr; PetscInt degrees[1000],ndegrees,npoints,two; PetscReal points[1000],weights[1000],interval[2]; PetscInt minpoints, maxpoints; PetscBool flg; PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Discretization tools test options",NULL);PetscCall(ierr); { ndegrees = 1000; degrees[0] = 0; degrees[1] = 1; degrees[2] = 2; PetscCall(PetscOptionsIntArray("-degrees","list of degrees to evaluate","",degrees,&ndegrees,&flg)); if (!flg) ndegrees = 3; npoints = 1000; points[0] = 0.0; points[1] = -0.5; points[2] = 1.0; PetscCall(PetscOptionsRealArray("-points","list of points at which to evaluate","",points,&npoints,&flg)); if (!flg) npoints = 3; two = 2; interval[0] = -1.; interval[1] = 1.; PetscCall(PetscOptionsRealArray("-interval","interval on which to construct quadrature","",interval,&two,NULL)); minpoints = 1; PetscCall(PetscOptionsInt("-minpoints","minimum points for thorough Gauss-Jacobi quadrature tests","",minpoints,&minpoints,NULL)); maxpoints = 30; #if defined(PETSC_USE_REAL_SINGLE) maxpoints = 5; #elif defined(PETSC_USE_REAL___FLOAT128) maxpoints = 20; /* just to make test faster */ #endif PetscCall(PetscOptionsInt("-maxpoints","maximum points for thorough Gauss-Jacobi quadrature tests","",maxpoints,&maxpoints,NULL)); } ierr = PetscOptionsEnd();PetscCall(ierr); PetscCall(CheckPoints("User-provided points",npoints,points,ndegrees,degrees)); PetscCall(PetscDTGaussQuadrature(npoints,interval[0],interval[1],points,weights)); PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Quadrature weights\n")); PetscCall(PetscRealView(npoints,weights,PETSC_VIEWER_STDOUT_WORLD)); { PetscReal a = interval[0],b = interval[1],zeroth,first,second; PetscInt i; zeroth = b - a; first = (b*b - a*a)/2; second = (b*b*b - a*a*a)/3; for (i=0; i= 2) { PetscCall(CheckJacobiQuadrature(i, 0., 0., PetscDTGaussLobattoJacobiQuadrature, 2*i-3)); PetscCall(CheckJacobiQuadrature(i, a1, b1, PetscDTGaussLobattoJacobiQuadrature, 2*i-3)); PetscCall(CheckJacobiQuadrature(i, a2, b2, PetscDTGaussLobattoJacobiQuadrature, 2*i-3)); } } } PetscCall(PetscFinalize()); return 0; } /*TEST test: suffix: 1 args: -degrees 1,2,3,4,5 -points 0,.2,-.5,.8,.9,1 -interval -.5,1 TEST*/