xref: /petsc/src/sys/tests/ex34.c (revision 74ac20f2b439a497e43401395505f0234d4acc29)
1c4762a1bSJed Brown static char help[] = "Tests IsInf/IsNan routines.\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,"%-32s -> %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 neg_zero = PetscRealConstant(-0.0);
18c4762a1bSJed Brown   PetscReal pos_zero = PetscRealConstant(+0.0);
19c4762a1bSJed Brown   PetscReal neg_one  = PetscRealConstant(-1.0);
20c4762a1bSJed Brown   PetscReal pos_one  = PetscRealConstant(+1.0);
21c4762a1bSJed Brown   PetscReal neg_inf  = neg_one/zero; /* -inf */
22c4762a1bSJed Brown   PetscReal pos_inf  = pos_one/zero; /* +inf */
23*74ac20f2SStefano Zampini   PetscReal x_nan    = zero2/zero;   /*  NaN */ /* some compilers may optimize out zero/zero and set x_nan = 1! */
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   PetscErrorCode ierr;
26c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   CALL(PetscIsInfReal(neg_zero));
29c4762a1bSJed Brown   CALL(PetscIsInfReal(pos_zero));
30c4762a1bSJed Brown   CALL(PetscIsInfReal(neg_one));
31c4762a1bSJed Brown   CALL(PetscIsInfReal(pos_one));
32c4762a1bSJed Brown   CALL(PetscIsInfReal(neg_inf));
33c4762a1bSJed Brown   CALL(PetscIsInfReal(pos_inf));
34c4762a1bSJed Brown   CALL(PetscIsInfReal(x_nan));
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   CALL(PetscIsNanReal(neg_zero));
37c4762a1bSJed Brown   CALL(PetscIsNanReal(pos_zero));
38c4762a1bSJed Brown   CALL(PetscIsNanReal(neg_one));
39c4762a1bSJed Brown   CALL(PetscIsNanReal(pos_one));
40c4762a1bSJed Brown   CALL(PetscIsNanReal(neg_inf));
41c4762a1bSJed Brown   CALL(PetscIsNanReal(pos_inf));
42c4762a1bSJed Brown   CALL(PetscIsNanReal(x_nan));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   ierr = PetscFinalize();
45c4762a1bSJed Brown   return ierr;
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /*TEST
49c4762a1bSJed Brown 
50c4762a1bSJed Brown    test:
51c4762a1bSJed Brown       output_file: output/ex34.out
52c4762a1bSJed Brown 
53c4762a1bSJed Brown TEST*/
54