1c4762a1bSJed Brown static char help[] = "Tests 1D discretization tools.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdt.h> 4c4762a1bSJed Brown #include <petscviewer.h> 5c4762a1bSJed Brown #include <petsc/private/petscimpl.h> 6c4762a1bSJed Brown #include <petsc/private/dtimpl.h> 7c4762a1bSJed Brown 8c4762a1bSJed Brown static PetscErrorCode CheckPoints(const char *name,PetscInt npoints,const PetscReal *points,PetscInt ndegrees,const PetscInt *degrees) 9c4762a1bSJed Brown { 10c4762a1bSJed Brown PetscReal *B,*D,*D2; 11c4762a1bSJed Brown PetscInt i,j; 12c4762a1bSJed Brown 13c4762a1bSJed Brown PetscFunctionBegin; 149566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(npoints*ndegrees,&B,npoints*ndegrees,&D,npoints*ndegrees,&D2)); 159566063dSJacob Faibussowitsch PetscCall(PetscDTLegendreEval(npoints,points,ndegrees,degrees,B,D,D2)); 169566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s\n",name)); 17c4762a1bSJed Brown for (i=0; i<npoints; i++) { 18c4762a1bSJed Brown for (j=0; j<ndegrees; j++) { 19c4762a1bSJed Brown PetscReal b,d,d2; 20c4762a1bSJed Brown b = B[i*ndegrees+j]; 21c4762a1bSJed Brown d = D[i*ndegrees+j]; 22c4762a1bSJed Brown d2 = D2[i*ndegrees+j]; 23c4762a1bSJed Brown if (PetscAbsReal(b) < PETSC_SMALL) b = 0; 24c4762a1bSJed Brown if (PetscAbsReal(d) < PETSC_SMALL) d = 0; 25c4762a1bSJed Brown if (PetscAbsReal(d2) < PETSC_SMALL) d2 = 0; 2663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"degree %" PetscInt_FMT " at %12.4g: B=%12.4g D=%12.4g D2=%12.4g\n",degrees[j],(double)points[i],(double)b,(double)d,(double)d2)); 27c4762a1bSJed Brown } 28c4762a1bSJed Brown } 299566063dSJacob Faibussowitsch PetscCall(PetscFree3(B,D,D2)); 30c4762a1bSJed Brown PetscFunctionReturn(0); 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33c4762a1bSJed Brown typedef PetscErrorCode(*quadratureFunc)(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal[],PetscReal[]); 34c4762a1bSJed Brown 35c4762a1bSJed Brown static PetscErrorCode CheckQuadrature_Basics(PetscInt npoints, PetscReal alpha, PetscReal beta, const PetscReal x[], const PetscReal w[]) 36c4762a1bSJed Brown { 37c4762a1bSJed Brown PetscInt i; 38c4762a1bSJed Brown 39c4762a1bSJed Brown PetscFunctionBegin; 40c4762a1bSJed Brown for (i = 1; i < npoints; i++) { 4163a3b9bcSJacob Faibussowitsch PetscCheck(x[i] > x[i-1],PETSC_COMM_SELF,PETSC_ERR_PLIB,"Quadrature points not monotonically increasing, %" PetscInt_FMT " points, alpha = %g, beta = %g, i = %" PetscInt_FMT ", x[i] = %g, x[i-1] = %g",npoints, (double) alpha, (double) beta, i, (double)x[i], (double)x[i-1]); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown for (i = 0; i < npoints; i++) { 4463a3b9bcSJacob Faibussowitsch PetscCheck(w[i] > 0.,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Quadrature weight not positive, %" PetscInt_FMT " points, alpha = %g, beta = %g, i = %" PetscInt_FMT ", w[i] = %g",npoints, (double) alpha, (double) beta, i, (double)w[i]); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown PetscFunctionReturn(0); 47c4762a1bSJed Brown } 48c4762a1bSJed Brown 49c4762a1bSJed Brown static PetscErrorCode CheckQuadrature(PetscInt npoints, PetscReal alpha, PetscReal beta, const PetscReal x[], const PetscReal w[], PetscInt nexact) 50c4762a1bSJed Brown { 51c4762a1bSJed Brown PetscInt i, j, k; 52c4762a1bSJed Brown PetscReal *Pi, *Pj; 53c4762a1bSJed Brown PetscReal eps; 54c4762a1bSJed Brown 55c4762a1bSJed Brown PetscFunctionBegin; 56c4762a1bSJed Brown eps = PETSC_SMALL; 579566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(npoints, &Pi, npoints, &Pj)); 58c4762a1bSJed Brown for (i = 0; i <= nexact; i++) { 599566063dSJacob Faibussowitsch PetscCall(PetscDTJacobiEval(npoints, alpha, beta, x, 1, &i, Pi, NULL, NULL)); 60c4762a1bSJed Brown for (j = i; j <= nexact - i; j++) { 61c4762a1bSJed Brown PetscReal I_quad = 0.; 62c4762a1bSJed Brown PetscReal I_exact = 0.; 63c4762a1bSJed Brown PetscReal err, tol; 649566063dSJacob Faibussowitsch PetscCall(PetscDTJacobiEval(npoints, alpha, beta, x, 1, &j, Pj, NULL, NULL)); 65c4762a1bSJed Brown 66c4762a1bSJed Brown tol = eps; 67c4762a1bSJed Brown if (i == j) { 68fbdc3dfeSToby Isaac PetscReal norm, norm2diff; 69fbdc3dfeSToby Isaac 70c4762a1bSJed Brown I_exact = PetscPowReal(2.0, alpha + beta + 1.) / (2.*i + alpha + beta + 1.); 71c4762a1bSJed Brown #if defined(PETSC_HAVE_LGAMMA) 72c4762a1bSJed Brown I_exact *= PetscExpReal(PetscLGamma(i + alpha + 1.) + PetscLGamma(i + beta + 1.) - (PetscLGamma(i + alpha + beta + 1.) + PetscLGamma(i + 1.))); 73c4762a1bSJed Brown #else 74c4762a1bSJed Brown { 75c4762a1bSJed Brown PetscInt ibeta = (PetscInt) beta; 76c4762a1bSJed Brown 7708401ef6SPierre Jolivet PetscCheck((PetscReal) ibeta == beta,PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable."); 78c4762a1bSJed Brown for (k = 0; k < ibeta; k++) I_exact *= (i + 1. + k) / (i + alpha + 1. + k); 79c4762a1bSJed Brown } 80c4762a1bSJed Brown #endif 81fbdc3dfeSToby Isaac 829566063dSJacob Faibussowitsch PetscCall(PetscDTJacobiNorm(alpha, beta, i, &norm)); 83fbdc3dfeSToby Isaac norm2diff = PetscAbsReal(norm*norm - I_exact); 84*1dca8a05SBarry Smith PetscCheck(norm2diff <= eps * I_exact,PETSC_COMM_SELF,PETSC_ERR_PLIB, "Jacobi norm error %g", (double) norm2diff); 85fbdc3dfeSToby Isaac 86c4762a1bSJed Brown tol = eps * I_exact; 87c4762a1bSJed Brown } 88c4762a1bSJed Brown for (k = 0; k < npoints; k++) I_quad += w[k] * (Pi[k] * Pj[k]); 89c4762a1bSJed Brown err = PetscAbsReal(I_exact - I_quad); 9063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"npoints %" PetscInt_FMT ", alpha %g, beta %g, i %" PetscInt_FMT ", j %" PetscInt_FMT ", exact %g, err %g\n", npoints, (double) alpha, (double) beta, i, j, (double) I_exact, (double) err)); 9163a3b9bcSJacob Faibussowitsch PetscCheck(err <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrectly integrated P_%" PetscInt_FMT " * P_%" PetscInt_FMT " using %" PetscInt_FMT " point rule with alpha = %g, beta = %g: exact %g, err %g", i, j, npoints, (double) alpha, (double) beta, (double) I_exact, (double) err); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown } 949566063dSJacob Faibussowitsch PetscCall(PetscFree2(Pi, Pj)); 95c4762a1bSJed Brown PetscFunctionReturn(0); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98c4762a1bSJed Brown static PetscErrorCode CheckJacobiQuadrature(PetscInt npoints, PetscReal alpha, PetscReal beta, quadratureFunc func, PetscInt nexact) 99c4762a1bSJed Brown { 100c4762a1bSJed Brown PetscReal *x, *w; 101c4762a1bSJed Brown 102c4762a1bSJed Brown PetscFunctionBegin; 1039566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(npoints, &x, npoints, &w)); 1049566063dSJacob Faibussowitsch PetscCall((*func)(npoints, -1., 1., alpha, beta, x, w)); 1059566063dSJacob Faibussowitsch PetscCall(CheckQuadrature_Basics(npoints, alpha, beta, x, w)); 1069566063dSJacob Faibussowitsch PetscCall(CheckQuadrature(npoints, alpha, beta, x, w, nexact)); 107c4762a1bSJed Brown #if defined(PETSCDTGAUSSIANQUADRATURE_EIG) 108c4762a1bSJed Brown /* compare methods of computing quadrature */ 109c4762a1bSJed Brown PetscDTGaussQuadratureNewton_Internal = (PetscBool) !PetscDTGaussQuadratureNewton_Internal; 110c4762a1bSJed Brown { 111c4762a1bSJed Brown PetscReal *x2, *w2; 112c4762a1bSJed Brown PetscReal eps; 113c4762a1bSJed Brown PetscInt i; 114c4762a1bSJed Brown 115c4762a1bSJed Brown eps = PETSC_SMALL; 1169566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(npoints, &x2, npoints, &w2)); 1179566063dSJacob Faibussowitsch PetscCall((*func)(npoints, -1., 1., alpha, beta, x2, w2)); 1189566063dSJacob Faibussowitsch PetscCall(CheckQuadrature_Basics(npoints, alpha, beta, x2, w2)); 1199566063dSJacob Faibussowitsch PetscCall(CheckQuadrature(npoints, alpha, beta, x2, w2, nexact)); 120c4762a1bSJed Brown for (i = 0; i < npoints; i++) { 121c4762a1bSJed Brown PetscReal xdiff, xtol, wdiff, wtol; 122c4762a1bSJed Brown 123c4762a1bSJed Brown xdiff = PetscAbsReal(x[i] - x2[i]); 124c4762a1bSJed Brown wdiff = PetscAbsReal(w[i] - w2[i]); 125c4762a1bSJed Brown xtol = eps * (1. + PetscMin(PetscAbsReal(x[i]),1. - PetscAbsReal(x[i]))); 126c4762a1bSJed Brown wtol = eps * (1. + w[i]); 12763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"npoints %" PetscInt_FMT ", alpha %g, beta %g, i %" PetscInt_FMT ", xdiff/xtol %g, wdiff/wtol %g\n", npoints, (double) alpha, (double) beta, i, (double)(xdiff/xtol), (double)(wdiff/wtol))); 12863a3b9bcSJacob Faibussowitsch PetscCheck(xdiff <= xtol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Mismatch quadrature point: %" PetscInt_FMT " points, alpha = %g, beta = %g, i = %" PetscInt_FMT ", xdiff = %g", npoints, (double) alpha, (double) beta, i, (double) xdiff); 12963a3b9bcSJacob Faibussowitsch PetscCheck(wdiff <= wtol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Mismatch quadrature weight: %" PetscInt_FMT " points, alpha = %g, beta = %g, i = %" PetscInt_FMT ", wdiff = %g", npoints, (double) alpha, (double) beta, i, (double) wdiff); 130c4762a1bSJed Brown } 1319566063dSJacob Faibussowitsch PetscCall(PetscFree2(x2, w2)); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown /* restore */ 134c4762a1bSJed Brown PetscDTGaussQuadratureNewton_Internal = (PetscBool) !PetscDTGaussQuadratureNewton_Internal; 135c4762a1bSJed Brown #endif 1369566063dSJacob Faibussowitsch PetscCall(PetscFree2(x, w)); 137c4762a1bSJed Brown PetscFunctionReturn(0); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown 140c4762a1bSJed Brown int main(int argc,char **argv) 141c4762a1bSJed Brown { 142c4762a1bSJed Brown PetscInt degrees[1000],ndegrees,npoints,two; 143c4762a1bSJed Brown PetscReal points[1000],weights[1000],interval[2]; 144c4762a1bSJed Brown PetscInt minpoints, maxpoints; 145c4762a1bSJed Brown PetscBool flg; 146c4762a1bSJed Brown 1479566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 148d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Discretization tools test options",NULL); 149c4762a1bSJed Brown { 150c4762a1bSJed Brown ndegrees = 1000; 151c4762a1bSJed Brown degrees[0] = 0; 152c4762a1bSJed Brown degrees[1] = 1; 153c4762a1bSJed Brown degrees[2] = 2; 1549566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-degrees","list of degrees to evaluate","",degrees,&ndegrees,&flg)); 155c4762a1bSJed Brown 156c4762a1bSJed Brown if (!flg) ndegrees = 3; 157c4762a1bSJed Brown npoints = 1000; 158c4762a1bSJed Brown points[0] = 0.0; 159c4762a1bSJed Brown points[1] = -0.5; 160c4762a1bSJed Brown points[2] = 1.0; 1619566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-points","list of points at which to evaluate","",points,&npoints,&flg)); 162c4762a1bSJed Brown 163c4762a1bSJed Brown if (!flg) npoints = 3; 164c4762a1bSJed Brown two = 2; 165c4762a1bSJed Brown interval[0] = -1.; 166c4762a1bSJed Brown interval[1] = 1.; 1679566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-interval","interval on which to construct quadrature","",interval,&two,NULL)); 168c4762a1bSJed Brown 169c4762a1bSJed Brown minpoints = 1; 1709566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-minpoints","minimum points for thorough Gauss-Jacobi quadrature tests","",minpoints,&minpoints,NULL)); 171c4762a1bSJed Brown maxpoints = 30; 172c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) 173c4762a1bSJed Brown maxpoints = 5; 174c4762a1bSJed Brown #elif defined(PETSC_USE_REAL___FLOAT128) 175c4762a1bSJed Brown maxpoints = 20; /* just to make test faster */ 176c4762a1bSJed Brown #endif 1779566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-maxpoints","maximum points for thorough Gauss-Jacobi quadrature tests","",maxpoints,&maxpoints,NULL)); 178c4762a1bSJed Brown } 179d0609cedSBarry Smith PetscOptionsEnd(); 1809566063dSJacob Faibussowitsch PetscCall(CheckPoints("User-provided points",npoints,points,ndegrees,degrees)); 181c4762a1bSJed Brown 1829566063dSJacob Faibussowitsch PetscCall(PetscDTGaussQuadrature(npoints,interval[0],interval[1],points,weights)); 1839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Quadrature weights\n")); 1849566063dSJacob Faibussowitsch PetscCall(PetscRealView(npoints,weights,PETSC_VIEWER_STDOUT_WORLD)); 185c4762a1bSJed Brown { 186c4762a1bSJed Brown PetscReal a = interval[0],b = interval[1],zeroth,first,second; 187c4762a1bSJed Brown PetscInt i; 188c4762a1bSJed Brown zeroth = b - a; 189c4762a1bSJed Brown first = (b*b - a*a)/2; 190c4762a1bSJed Brown second = (b*b*b - a*a*a)/3; 191c4762a1bSJed Brown for (i=0; i<npoints; i++) { 192c4762a1bSJed Brown zeroth -= weights[i]; 193c4762a1bSJed Brown first -= weights[i] * points[i]; 194c4762a1bSJed Brown second -= weights[i] * PetscSqr(points[i]); 195c4762a1bSJed Brown } 196c4762a1bSJed Brown if (PetscAbs(zeroth) < 1e-10) zeroth = 0.; 197c4762a1bSJed Brown if (PetscAbs(first) < 1e-10) first = 0.; 198c4762a1bSJed Brown if (PetscAbs(second) < 1e-10) second = 0.; 1999566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Moment error: zeroth=%g, first=%g, second=%g\n",(double)(-zeroth),(double)(-first),(double)(-second))); 200c4762a1bSJed Brown } 2019566063dSJacob Faibussowitsch PetscCall(CheckPoints("Gauss points",npoints,points,ndegrees,degrees)); 202c4762a1bSJed Brown { 203c4762a1bSJed Brown PetscInt i; 204c4762a1bSJed Brown 205c4762a1bSJed Brown for (i = minpoints; i <= maxpoints; i++) { 206c4762a1bSJed Brown PetscReal a1, b1, a2, b2; 207c4762a1bSJed Brown 208c4762a1bSJed Brown #if defined(PETSC_HAVE_LGAMMA) 209c4762a1bSJed Brown a1 = -0.6; 210c4762a1bSJed Brown b1 = 1.1; 211c4762a1bSJed Brown a2 = 2.2; 212c4762a1bSJed Brown b2 = -0.6; 213c4762a1bSJed Brown #else 214c4762a1bSJed Brown a1 = 0.; 215c4762a1bSJed Brown b1 = 1.; 216c4762a1bSJed Brown a2 = 2.; 217c4762a1bSJed Brown b2 = 0.; 218c4762a1bSJed Brown #endif 2199566063dSJacob Faibussowitsch PetscCall(CheckJacobiQuadrature(i, 0., 0., PetscDTGaussJacobiQuadrature, 2*i-1)); 2209566063dSJacob Faibussowitsch PetscCall(CheckJacobiQuadrature(i, a1, b1, PetscDTGaussJacobiQuadrature, 2*i-1)); 2219566063dSJacob Faibussowitsch PetscCall(CheckJacobiQuadrature(i, a2, b2, PetscDTGaussJacobiQuadrature, 2*i-1)); 222c4762a1bSJed Brown if (i >= 2) { 2239566063dSJacob Faibussowitsch PetscCall(CheckJacobiQuadrature(i, 0., 0., PetscDTGaussLobattoJacobiQuadrature, 2*i-3)); 2249566063dSJacob Faibussowitsch PetscCall(CheckJacobiQuadrature(i, a1, b1, PetscDTGaussLobattoJacobiQuadrature, 2*i-3)); 2259566063dSJacob Faibussowitsch PetscCall(CheckJacobiQuadrature(i, a2, b2, PetscDTGaussLobattoJacobiQuadrature, 2*i-3)); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown } 228c4762a1bSJed Brown } 2299566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 230b122ec5aSJacob Faibussowitsch return 0; 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233c4762a1bSJed Brown /*TEST 234c4762a1bSJed Brown test: 235c4762a1bSJed Brown suffix: 1 236c4762a1bSJed Brown args: -degrees 1,2,3,4,5 -points 0,.2,-.5,.8,.9,1 -interval -.5,1 237c4762a1bSJed Brown TEST*/ 238