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