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