1*c4762a1bSJed Brown static char help[] = "Tests quadrature.\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscdt.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown static void func1(PetscReal x, PetscReal *val) 6*c4762a1bSJed Brown { 7*c4762a1bSJed Brown *val = x*PetscLogReal(1+x); 8*c4762a1bSJed Brown } 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown static void func2(PetscReal x, PetscReal *val) 11*c4762a1bSJed Brown { 12*c4762a1bSJed Brown *val = x*x*PetscAtanReal(x); 13*c4762a1bSJed Brown } 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown static void func3(PetscReal x, PetscReal *val) 16*c4762a1bSJed Brown { 17*c4762a1bSJed Brown *val = PetscExpReal(x)*PetscCosReal(x); 18*c4762a1bSJed Brown } 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown static void func4(PetscReal x, PetscReal *val) 21*c4762a1bSJed Brown { 22*c4762a1bSJed Brown const PetscReal u = PetscSqrtReal(2.0 + x*x); 23*c4762a1bSJed Brown *val = PetscAtanReal(u)/((1.0 + x*x)*u); 24*c4762a1bSJed Brown } 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown static void func5(PetscReal x, PetscReal *val) 27*c4762a1bSJed Brown { 28*c4762a1bSJed Brown if (x == 0.0) *val = 0.0; 29*c4762a1bSJed Brown else *val = PetscSqrtReal(x)*PetscLogReal(x); 30*c4762a1bSJed Brown } 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown static void func6(PetscReal x, PetscReal *val) 33*c4762a1bSJed Brown { 34*c4762a1bSJed Brown *val = PetscSqrtReal(1-x*x); 35*c4762a1bSJed Brown } 36*c4762a1bSJed Brown 37*c4762a1bSJed Brown static void func7(PetscReal x, PetscReal *val) 38*c4762a1bSJed Brown { 39*c4762a1bSJed Brown if (x == 1.0) *val = PETSC_INFINITY; 40*c4762a1bSJed Brown else *val = PetscSqrtReal(x)/PetscSqrtReal(1-x*x); 41*c4762a1bSJed Brown } 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown static void func8(PetscReal x, PetscReal *val) 44*c4762a1bSJed Brown { 45*c4762a1bSJed Brown if (x == 0.0) *val = PETSC_INFINITY; 46*c4762a1bSJed Brown else *val = PetscLogReal(x)*PetscLogReal(x); 47*c4762a1bSJed Brown } 48*c4762a1bSJed Brown 49*c4762a1bSJed Brown static void func9(PetscReal x, PetscReal *val) 50*c4762a1bSJed Brown { 51*c4762a1bSJed Brown *val = PetscLogReal(PetscCosReal(x)); 52*c4762a1bSJed Brown } 53*c4762a1bSJed Brown 54*c4762a1bSJed Brown static void func10(PetscReal x, PetscReal *val) 55*c4762a1bSJed Brown { 56*c4762a1bSJed Brown if (x == 0.0) *val = 0.0; 57*c4762a1bSJed Brown else if (x == 1.0) *val = PETSC_INFINITY; 58*c4762a1bSJed Brown *val = PetscSqrtReal(PetscTanReal(x)); 59*c4762a1bSJed Brown } 60*c4762a1bSJed Brown 61*c4762a1bSJed Brown static void func11(PetscReal x, PetscReal *val) 62*c4762a1bSJed Brown { 63*c4762a1bSJed Brown *val = 1/(1-2*x+2*x*x); 64*c4762a1bSJed Brown } 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown static void func12(PetscReal x, PetscReal *val) 67*c4762a1bSJed Brown { 68*c4762a1bSJed Brown if (x == 0.0) *val = 0.0; 69*c4762a1bSJed Brown else if (x == 1.0) *val = PETSC_INFINITY; 70*c4762a1bSJed Brown else *val = PetscExpReal(1-1/x)/PetscSqrtReal(x*x*x-x*x*x*x); 71*c4762a1bSJed Brown } 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown static void func13(PetscReal x, PetscReal *val) 74*c4762a1bSJed Brown { 75*c4762a1bSJed Brown if (x == 0.0) *val = 0.0; 76*c4762a1bSJed Brown else if (x == 1.0) *val = 1.0; 77*c4762a1bSJed Brown else *val = PetscExpReal(-(1/x-1)*(1/x-1)/2)/(x*x); 78*c4762a1bSJed Brown } 79*c4762a1bSJed Brown 80*c4762a1bSJed Brown static void func14(PetscReal x, PetscReal *val) 81*c4762a1bSJed Brown { 82*c4762a1bSJed Brown if (x == 0.0) *val = 0.0; 83*c4762a1bSJed Brown else if (x == 1.0) *val = 1.0; 84*c4762a1bSJed Brown else *val = PetscExpReal(1-1/x)*PetscCosReal(1/x-1)/(x*x); 85*c4762a1bSJed Brown } 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown int main(int argc, char **argv) 88*c4762a1bSJed Brown { 89*c4762a1bSJed Brown #if PETSC_SCALAR_SIZE == 32 90*c4762a1bSJed Brown PetscInt digits = 7; 91*c4762a1bSJed Brown #elif PETSC_SCALAR_SIZE == 64 92*c4762a1bSJed Brown PetscInt digits = 14; 93*c4762a1bSJed Brown #else 94*c4762a1bSJed Brown PetscInt digits = 14; 95*c4762a1bSJed Brown #endif 96*c4762a1bSJed Brown /* for some reason in __float128 precision it cannot get more accuracy for some of the integrals */ 97*c4762a1bSJed Brown #if defined(PETSC_USE_REAL___FLOAT128) 98*c4762a1bSJed Brown const PetscReal epsilon = 2.2204460492503131e-16; 99*c4762a1bSJed Brown #else 100*c4762a1bSJed Brown const PetscReal epsilon = 2500.*PETSC_MACHINE_EPSILON; 101*c4762a1bSJed Brown #endif 102*c4762a1bSJed Brown const PetscReal bounds[28] = 103*c4762a1bSJed Brown { 104*c4762a1bSJed Brown 0.0, 1.0, 105*c4762a1bSJed Brown 0.0, 1.0, 106*c4762a1bSJed Brown 0.0, PETSC_PI/2., 107*c4762a1bSJed Brown 0.0, 1.0, 108*c4762a1bSJed Brown 0.0, 1.0, 109*c4762a1bSJed Brown 0.0, 1.0, 110*c4762a1bSJed Brown 0.0, 1.0, 111*c4762a1bSJed Brown 0.0, 1.0, 112*c4762a1bSJed Brown 0.0, PETSC_PI/2., 113*c4762a1bSJed Brown 0.0, PETSC_PI/2., 114*c4762a1bSJed Brown 0.0, 1.0, 115*c4762a1bSJed Brown 0.0, 1.0, 116*c4762a1bSJed Brown 0.0, 1.0, 117*c4762a1bSJed Brown 0.0, 1.0 118*c4762a1bSJed Brown }; 119*c4762a1bSJed Brown const PetscReal analytic[14] = 120*c4762a1bSJed Brown { 121*c4762a1bSJed Brown 0.250000000000000, 122*c4762a1bSJed Brown 0.210657251225806988108092302182988001695680805674, 123*c4762a1bSJed Brown 1.905238690482675827736517833351916563195085437332, 124*c4762a1bSJed Brown 0.514041895890070761397629739576882871630921844127, 125*c4762a1bSJed Brown -.444444444444444444444444444444444444444444444444, 126*c4762a1bSJed Brown 0.785398163397448309615660845819875721049292349843, 127*c4762a1bSJed Brown 1.198140234735592207439922492280323878227212663216, 128*c4762a1bSJed Brown 2.000000000000000000000000000000000000000000000000, 129*c4762a1bSJed Brown -1.08879304515180106525034444911880697366929185018, 130*c4762a1bSJed Brown 2.221441469079183123507940495030346849307310844687, 131*c4762a1bSJed Brown 1.570796326794896619231321691639751442098584699687, 132*c4762a1bSJed Brown 1.772453850905516027298167483341145182797549456122, 133*c4762a1bSJed Brown 1.253314137315500251207882642405522626503493370304, 134*c4762a1bSJed Brown 0.500000000000000000000000000000000000000000000000 135*c4762a1bSJed Brown }; 136*c4762a1bSJed Brown void (*funcs[14])(PetscReal, PetscReal *) = {func1, func2, func3, func4, func5, func6, func7, func8, func9, func10, func11, func12, func13, func14}; 137*c4762a1bSJed Brown PetscInt f; 138*c4762a1bSJed Brown PetscErrorCode ierr; 139*c4762a1bSJed Brown 140*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 141*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Test Options","none");CHKERRQ(ierr); 142*c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-digits", "The number of significant digits for the integral","ex3.c",digits,&digits,NULL,1);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = PetscOptionsEnd(); 144*c4762a1bSJed Brown 145*c4762a1bSJed Brown /* Integrate each function */ 146*c4762a1bSJed Brown for (f = 0; f < 14; ++f) { 147*c4762a1bSJed Brown PetscReal integral; 148*c4762a1bSJed Brown 149*c4762a1bSJed Brown /* These can only be integrated accuractely using MPFR */ 150*c4762a1bSJed Brown if ((f == 6) || (f == 7) || (f == 9) || (f == 11)) continue; 151*c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) 152*c4762a1bSJed Brown if (f == 8) continue; 153*c4762a1bSJed Brown #endif 154*c4762a1bSJed Brown ierr = PetscDTTanhSinhIntegrate(funcs[f], bounds[f*2+0], bounds[f*2+1], digits, &integral);CHKERRQ(ierr); 155*c4762a1bSJed Brown if (PetscAbsReal(integral - analytic[f]) > PetscMax(epsilon, PetscPowRealInt(10.0, -digits)) || PetscIsInfOrNanScalar(integral - analytic[f])) { 156*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "The integral of func%2d is wrong: %g (%g)\n", f+1, (double)integral, (double) PetscAbsReal(integral - analytic[f]));CHKERRQ(ierr); 157*c4762a1bSJed Brown } 158*c4762a1bSJed Brown } 159*c4762a1bSJed Brown #if defined(PETSC_HAVE_MPFR) 160*c4762a1bSJed Brown for (f = 0; f < 14; ++f) { 161*c4762a1bSJed Brown PetscReal integral; 162*c4762a1bSJed Brown 163*c4762a1bSJed Brown ierr = PetscDTTanhSinhIntegrateMPFR(funcs[f], bounds[f*2+0], bounds[f*2+1], digits, &integral);CHKERRQ(ierr); 164*c4762a1bSJed Brown if (PetscAbsReal(integral - analytic[f]) > PetscPowRealInt(10.0, -digits)) { 165*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "The integral of func%2d is wrong: %g (%g)\n", f+1, (double)integral, (double)PetscAbsReal(integral - analytic[f]));CHKERRQ(ierr); 166*c4762a1bSJed Brown } 167*c4762a1bSJed Brown } 168*c4762a1bSJed Brown #endif 169*c4762a1bSJed Brown ierr = PetscFinalize(); 170*c4762a1bSJed Brown return ierr; 171*c4762a1bSJed Brown } 172*c4762a1bSJed Brown 173*c4762a1bSJed Brown /*TEST 174*c4762a1bSJed Brown test: 175*c4762a1bSJed Brown suffix: 0 176*c4762a1bSJed Brown TEST*/ 177