1c4762a1bSJed Brown static char help[] = "Tests quadrature.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdt.h> 4c4762a1bSJed Brown 59371c9d4SSatish Balay static void func1(const PetscReal a[], void *dummy, PetscReal *val) { 6d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 7c4762a1bSJed Brown *val = x * PetscLogReal(1 + x); 8c4762a1bSJed Brown } 9c4762a1bSJed Brown 109371c9d4SSatish Balay static void func2(const PetscReal a[], void *dummy, PetscReal *val) { 11d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 12c4762a1bSJed Brown *val = x * x * PetscAtanReal(x); 13c4762a1bSJed Brown } 14c4762a1bSJed Brown 159371c9d4SSatish Balay static void func3(const PetscReal a[], void *dummy, PetscReal *val) { 16d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 17c4762a1bSJed Brown *val = PetscExpReal(x) * PetscCosReal(x); 18c4762a1bSJed Brown } 19c4762a1bSJed Brown 209371c9d4SSatish Balay static void func4(const PetscReal a[], void *dummy, PetscReal *val) { 21d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 22c4762a1bSJed Brown const PetscReal u = PetscSqrtReal(2.0 + x * x); 23c4762a1bSJed Brown *val = PetscAtanReal(u) / ((1.0 + x * x) * u); 24c4762a1bSJed Brown } 25c4762a1bSJed Brown 269371c9d4SSatish Balay static void func5(const PetscReal a[], void *dummy, PetscReal *val) { 27d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 28c4762a1bSJed Brown if (x == 0.0) *val = 0.0; 29c4762a1bSJed Brown else *val = PetscSqrtReal(x) * PetscLogReal(x); 30c4762a1bSJed Brown } 31c4762a1bSJed Brown 329371c9d4SSatish Balay static void func6(const PetscReal a[], void *dummy, PetscReal *val) { 33d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 34c4762a1bSJed Brown *val = PetscSqrtReal(1 - x * x); 35c4762a1bSJed Brown } 36c4762a1bSJed Brown 379371c9d4SSatish Balay static void func7(const PetscReal a[], void *dummy, PetscReal *val) { 38d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 39c4762a1bSJed Brown if (x == 1.0) *val = PETSC_INFINITY; 40c4762a1bSJed Brown else *val = PetscSqrtReal(x) / PetscSqrtReal(1 - x * x); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 439371c9d4SSatish Balay static void func8(const PetscReal a[], void *dummy, PetscReal *val) { 44d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 45c4762a1bSJed Brown if (x == 0.0) *val = PETSC_INFINITY; 46c4762a1bSJed Brown else *val = PetscLogReal(x) * PetscLogReal(x); 47c4762a1bSJed Brown } 48c4762a1bSJed Brown 499371c9d4SSatish Balay static void func9(const PetscReal x[], void *dummy, PetscReal *val) { 50d6685f55SMatthew G. Knepley *val = PetscLogReal(PetscCosReal(x[0])); 51c4762a1bSJed Brown } 52c4762a1bSJed Brown 539371c9d4SSatish Balay static void func10(const PetscReal a[], void *dummy, PetscReal *val) { 54d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 55c4762a1bSJed Brown if (x == 0.0) *val = 0.0; 56c4762a1bSJed Brown else if (x == 1.0) *val = PETSC_INFINITY; 57c4762a1bSJed Brown *val = PetscSqrtReal(PetscTanReal(x)); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 609371c9d4SSatish Balay static void func11(const PetscReal a[], void *dummy, PetscReal *val) { 61d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 62c4762a1bSJed Brown *val = 1 / (1 - 2 * x + 2 * x * x); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 659371c9d4SSatish Balay static void func12(const PetscReal a[], void *dummy, PetscReal *val) { 66d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 67c4762a1bSJed Brown if (x == 0.0) *val = 0.0; 68c4762a1bSJed Brown else if (x == 1.0) *val = PETSC_INFINITY; 69c4762a1bSJed Brown else *val = PetscExpReal(1 - 1 / x) / PetscSqrtReal(x * x * x - x * x * x * x); 70c4762a1bSJed Brown } 71c4762a1bSJed Brown 729371c9d4SSatish Balay static void func13(const PetscReal a[], void *dummy, PetscReal *val) { 73d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 74c4762a1bSJed Brown if (x == 0.0) *val = 0.0; 75c4762a1bSJed Brown else if (x == 1.0) *val = 1.0; 76c4762a1bSJed Brown else *val = PetscExpReal(-(1 / x - 1) * (1 / x - 1) / 2) / (x * x); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 799371c9d4SSatish Balay static void func14(const PetscReal a[], void *dummy, PetscReal *val) { 80d6685f55SMatthew G. Knepley const PetscReal x = a[0]; 81c4762a1bSJed Brown if (x == 0.0) *val = 0.0; 82c4762a1bSJed Brown else if (x == 1.0) *val = 1.0; 83c4762a1bSJed Brown else *val = PetscExpReal(1 - 1 / x) * PetscCosReal(1 / x - 1) / (x * x); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 869371c9d4SSatish Balay int main(int argc, char **argv) { 87c4762a1bSJed Brown #if PETSC_SCALAR_SIZE == 32 88c4762a1bSJed Brown PetscInt digits = 7; 89c4762a1bSJed Brown #elif PETSC_SCALAR_SIZE == 64 90c4762a1bSJed Brown PetscInt digits = 14; 91c4762a1bSJed Brown #else 92c4762a1bSJed Brown PetscInt digits = 14; 93c4762a1bSJed Brown #endif 94c4762a1bSJed Brown /* for some reason in __float128 precision it cannot get more accuracy for some of the integrals */ 95c4762a1bSJed Brown #if defined(PETSC_USE_REAL___FLOAT128) 96c4762a1bSJed Brown const PetscReal epsilon = 2.2204460492503131e-16; 97c4762a1bSJed Brown #else 98c4762a1bSJed Brown const PetscReal epsilon = 2500. * PETSC_MACHINE_EPSILON; 99c4762a1bSJed Brown #endif 1009371c9d4SSatish Balay const PetscReal bounds[28] = {0.0, 1.0, 0.0, 1.0, 0.0, PETSC_PI / 2., 0.0, 1.0, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0, 0.0, PETSC_PI / 2., 0.0, PETSC_PI / 2., 0.0, 1.0, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0}; 1019371c9d4SSatish Balay const PetscReal analytic[14] = {0.250000000000000, 0.210657251225806988108092302182988001695680805674, 1.905238690482675827736517833351916563195085437332, 0.514041895890070761397629739576882871630921844127, -.444444444444444444444444444444444444444444444444, 0.785398163397448309615660845819875721049292349843, 1.198140234735592207439922492280323878227212663216, 2.000000000000000000000000000000000000000000000000, -1.08879304515180106525034444911880697366929185018, 2.221441469079183123507940495030346849307310844687, 1.570796326794896619231321691639751442098584699687, 1.772453850905516027298167483341145182797549456122, 1.253314137315500251207882642405522626503493370304, 0.500000000000000000000000000000000000000000000000}; 102d6685f55SMatthew G. Knepley void (*funcs[14])(const PetscReal[], void *, PetscReal *) = {func1, func2, func3, func4, func5, func6, func7, func8, func9, func10, func11, func12, func13, func14}; 103c4762a1bSJed Brown PetscInt f; 104c4762a1bSJed Brown 105327415f7SBarry Smith PetscFunctionBeginUser; 1069566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 107d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, "", "Test Options", "none"); 1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-digits", "The number of significant digits for the integral", "ex3.c", digits, &digits, NULL, 1)); 109d0609cedSBarry Smith PetscOptionsEnd(); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* Integrate each function */ 112c4762a1bSJed Brown for (f = 0; f < 14; ++f) { 113c4762a1bSJed Brown PetscReal integral; 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* These can only be integrated accuractely using MPFR */ 116c4762a1bSJed Brown if ((f == 6) || (f == 7) || (f == 9) || (f == 11)) continue; 117c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) 118c4762a1bSJed Brown if (f == 8) continue; 119c4762a1bSJed Brown #endif 1209566063dSJacob Faibussowitsch PetscCall(PetscDTTanhSinhIntegrate(funcs[f], bounds[f * 2 + 0], bounds[f * 2 + 1], digits, NULL, &integral)); 121c4762a1bSJed Brown if (PetscAbsReal(integral - analytic[f]) > PetscMax(epsilon, PetscPowRealInt(10.0, -digits)) || PetscIsInfOrNanScalar(integral - analytic[f])) { 12263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "The integral of func%2" PetscInt_FMT " is wrong: %g (%g)\n", f + 1, (double)integral, (double)PetscAbsReal(integral - analytic[f]))); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown } 125c4762a1bSJed Brown #if defined(PETSC_HAVE_MPFR) 126c4762a1bSJed Brown for (f = 0; f < 14; ++f) { 127c4762a1bSJed Brown PetscReal integral; 128c4762a1bSJed Brown 1299566063dSJacob Faibussowitsch PetscCall(PetscDTTanhSinhIntegrateMPFR(funcs[f], bounds[f * 2 + 0], bounds[f * 2 + 1], digits, NULL, &integral)); 130*48a46eb9SPierre Jolivet if (PetscAbsReal(integral - analytic[f]) > PetscPowRealInt(10.0, -digits)) PetscCall(PetscPrintf(PETSC_COMM_SELF, "The integral of func%2d is wrong: %g (%g)\n", f + 1, (double)integral, (double)PetscAbsReal(integral - analytic[f]))); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown #endif 1339566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 134b122ec5aSJacob Faibussowitsch return 0; 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown /*TEST 138c4762a1bSJed Brown test: 139c4762a1bSJed Brown suffix: 0 140c4762a1bSJed Brown TEST*/ 141