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