xref: /petsc/src/sys/tests/ex34.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
774ac20f2SStefano Zampini PETSC_INTERN PetscReal zero2;
874ac20f2SStefano Zampini PetscReal zero2 = 0;
9c4762a1bSJed Brown 
105f80ce2aSJacob Faibussowitsch #define CALL(call) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%-32s -> %s\n",#call,(call)?"True":"False"))
11c4762a1bSJed Brown 
12c4762a1bSJed Brown int main(int argc, char **argv) {
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   PetscReal neg_zero = PetscRealConstant(-0.0);
15c4762a1bSJed Brown   PetscReal pos_zero = PetscRealConstant(+0.0);
16c4762a1bSJed Brown   PetscReal neg_one  = PetscRealConstant(-1.0);
17c4762a1bSJed Brown   PetscReal pos_one  = PetscRealConstant(+1.0);
18c4762a1bSJed Brown   PetscReal neg_inf  = neg_one/zero; /* -inf */
19c4762a1bSJed Brown   PetscReal pos_inf  = pos_one/zero; /* +inf */
2074ac20f2SStefano Zampini   PetscReal x_nan    = zero2/zero;   /*  NaN */ /* some compilers may optimize out zero/zero and set x_nan = 1! */
21c4762a1bSJed Brown 
22*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   CALL(PetscIsInfReal(neg_zero));
25c4762a1bSJed Brown   CALL(PetscIsInfReal(pos_zero));
26c4762a1bSJed Brown   CALL(PetscIsInfReal(neg_one));
27c4762a1bSJed Brown   CALL(PetscIsInfReal(pos_one));
28c4762a1bSJed Brown   CALL(PetscIsInfReal(neg_inf));
29c4762a1bSJed Brown   CALL(PetscIsInfReal(pos_inf));
30c4762a1bSJed Brown   CALL(PetscIsInfReal(x_nan));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   CALL(PetscIsNanReal(neg_zero));
33c4762a1bSJed Brown   CALL(PetscIsNanReal(pos_zero));
34c4762a1bSJed Brown   CALL(PetscIsNanReal(neg_one));
35c4762a1bSJed Brown   CALL(PetscIsNanReal(pos_one));
36c4762a1bSJed Brown   CALL(PetscIsNanReal(neg_inf));
37c4762a1bSJed Brown   CALL(PetscIsNanReal(pos_inf));
38c4762a1bSJed Brown   CALL(PetscIsNanReal(x_nan));
39c4762a1bSJed Brown 
40*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
41*b122ec5aSJacob Faibussowitsch   return 0;
42c4762a1bSJed Brown }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown /*TEST
45c4762a1bSJed Brown 
46c4762a1bSJed Brown    test:
47c4762a1bSJed Brown       output_file: output/ex34.out
48c4762a1bSJed Brown 
49c4762a1bSJed Brown TEST*/
50