xref: /petsc/src/sys/tests/ex39.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 static char help[] = "Tests PetscIsCloseAtTol() routine.\n";
2 
3 #include <petscsys.h>
4 
5 PETSC_INTERN PetscReal zero;
6 PetscReal              zero = 0;
7 PETSC_INTERN PetscReal zero2;
8 PetscReal              zero2 = 0;
9 
10 #define CALL(call) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s -> %s\n", #call, (call) ? "True" : "False"))
11 
12 int main(int argc, char **argv) {
13   PetscReal eps      = PETSC_MACHINE_EPSILON;
14   PetscReal neg_zero = PetscRealConstant(-0.0);
15   PetscReal pos_zero = PetscRealConstant(+0.0);
16   PetscReal neg_one  = PetscRealConstant(-1.0);
17   PetscReal pos_one  = PetscRealConstant(+1.0);
18   PetscReal neg_inf  = neg_one / zero;          /* -inf */
19   PetscReal pos_inf  = pos_one / zero;          /* +inf */
20   PetscReal x_nan    = zero2 / zero; /*  NaN */ /* some compilers may optimize out zero/zero and set x_nan = 1! */
21 
22   PetscFunctionBeginUser;
23   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
24 
25   CALL(PetscIsCloseAtTol(pos_zero, neg_zero, 0, 0));
26   CALL(PetscIsCloseAtTol(pos_one, pos_one, 0, 0));
27   CALL(PetscIsCloseAtTol(pos_one, neg_one, 0, 0));
28   CALL(PetscIsCloseAtTol(pos_one, neg_one, 0, 2));
29   CALL(PetscIsCloseAtTol(pos_one, neg_one, 2, 0));
30 
31   CALL(PetscIsCloseAtTol(pos_one + eps, pos_one, 0, 0));
32   CALL(PetscIsCloseAtTol(pos_one - eps, pos_one, 0, 0));
33   CALL(PetscIsCloseAtTol(pos_one + eps, pos_one, 0, 0));
34   CALL(PetscIsCloseAtTol(pos_one - eps, pos_one, 0, 0));
35 
36   CALL(PetscIsCloseAtTol(pos_one + eps, pos_one, 0, eps));
37   CALL(PetscIsCloseAtTol(pos_one - eps, pos_one, 0, eps));
38   CALL(PetscIsCloseAtTol(pos_one + eps, pos_one, eps, 0));
39   CALL(PetscIsCloseAtTol(pos_one - eps, pos_one, eps, 0));
40 
41   CALL(PetscIsCloseAtTol(pos_one + 2 * eps, pos_one, eps, 0));
42   CALL(PetscIsCloseAtTol(pos_one - 2 * eps, pos_one, eps, 0));
43   CALL(PetscIsCloseAtTol(pos_one + 2 * eps, pos_one, 0, eps));
44   CALL(PetscIsCloseAtTol(pos_one - 2 * eps, pos_one, 0, eps));
45 
46   CALL(PetscIsCloseAtTol(neg_inf, neg_zero, 2, 2));
47   CALL(PetscIsCloseAtTol(neg_inf, pos_zero, 2, 2));
48   CALL(PetscIsCloseAtTol(neg_inf, neg_one, 2, 2));
49   CALL(PetscIsCloseAtTol(neg_inf, pos_one, 2, 2));
50   CALL(PetscIsCloseAtTol(neg_inf, neg_inf, 2, 2));
51   CALL(PetscIsCloseAtTol(neg_inf, pos_inf, 2, 2));
52   CALL(PetscIsCloseAtTol(neg_inf, x_nan, 2, 2));
53 
54   CALL(PetscIsCloseAtTol(pos_inf, neg_zero, 2, 2));
55   CALL(PetscIsCloseAtTol(pos_inf, pos_zero, 2, 2));
56   CALL(PetscIsCloseAtTol(pos_inf, neg_one, 2, 2));
57   CALL(PetscIsCloseAtTol(pos_inf, pos_one, 2, 2));
58   CALL(PetscIsCloseAtTol(pos_inf, neg_inf, 2, 2));
59   CALL(PetscIsCloseAtTol(pos_inf, pos_inf, 2, 2));
60   CALL(PetscIsCloseAtTol(pos_inf, x_nan, 2, 2));
61 
62   CALL(PetscIsCloseAtTol(x_nan, neg_zero, 2, 2));
63   CALL(PetscIsCloseAtTol(x_nan, pos_zero, 2, 2));
64   CALL(PetscIsCloseAtTol(x_nan, neg_one, 2, 2));
65   CALL(PetscIsCloseAtTol(x_nan, pos_one, 2, 2));
66   CALL(PetscIsCloseAtTol(x_nan, neg_inf, 2, 2));
67   CALL(PetscIsCloseAtTol(x_nan, pos_inf, 2, 2));
68   CALL(PetscIsCloseAtTol(x_nan, x_nan, 2, 2));
69 
70   PetscCall(PetscFinalize());
71   return 0;
72 }
73 
74 /*TEST
75 
76    test:
77       output_file: output/ex39.out
78 
79 TEST*/
80