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