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