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