10c567f5aSSatish Balay #define PETSC_SKIP_COMPLEX 2532fbc7fSSatish Balay #include <petscsys.h> 3532fbc7fSSatish Balay /*@C 4*811af0c4SBarry Smith PetscIsNormalReal - Returns `PETSC_TRUE` if the input value satisfies `isnormal()` 5532fbc7fSSatish Balay 6532fbc7fSSatish Balay Input Parameter: 78b49ba18SBarry Smith . a - the PetscReal Value 8532fbc7fSSatish Balay 9*811af0c4SBarry Smith Developer Notes: 10*811af0c4SBarry Smith Uses the C99 standard `isnormal()` on systems where they exist. 11*811af0c4SBarry Smith 12*811af0c4SBarry Smith Uses `isnormalq()` with `__float128` 13*811af0c4SBarry Smith 14*811af0c4SBarry Smith Otherwise always returns true 158b49ba18SBarry Smith 168b49ba18SBarry Smith Level: beginner 17*811af0c4SBarry Smith 18*811af0c4SBarry Smith .seealso: `PetscIsInfReal()`, `PetscIsNanReal()` 198b49ba18SBarry Smith @*/ 20570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 219371c9d4SSatish Balay PetscBool PetscIsNormalReal(PetscReal a) { 228fa295daSBarry Smith return PETSC_TRUE; 238b49ba18SBarry Smith } 248b49ba18SBarry Smith #elif defined(PETSC_HAVE_ISNORMAL) 259371c9d4SSatish Balay PetscBool PetscIsNormalReal(PetscReal a) { 26bb88209dSBarry Smith return isnormal(a) ? PETSC_TRUE : PETSC_FALSE; 278b49ba18SBarry Smith } 288b49ba18SBarry Smith #else 299371c9d4SSatish Balay PetscBool PetscIsNormalReal(PetscReal a) { 308b49ba18SBarry Smith return PETSC_TRUE; 318b49ba18SBarry Smith } 328b49ba18SBarry Smith #endif 338b49ba18SBarry Smith 348b49ba18SBarry Smith /*@C 359f4f8022SLisandro Dalcin PetscIsInfReal - Returns whether the input is an infinity value. 368b49ba18SBarry Smith 378b49ba18SBarry Smith Input Parameter: 388b49ba18SBarry Smith . a - the floating point number 39532fbc7fSSatish Balay 40*811af0c4SBarry Smith Developer Notes: 41*811af0c4SBarry Smith Uses the C99 standard `isinf()` on systems where it exists. 42*811af0c4SBarry Smith 43*811af0c4SBarry Smith Otherwise uses (a && a/2 == a), note that some optimizing compilers compile out this form, thus removing the check. 44532fbc7fSSatish Balay 45532fbc7fSSatish Balay Level: beginner 46*811af0c4SBarry Smith 47*811af0c4SBarry Smith .seealso: `PetscIsNormalReal()`, `PetscIsNanReal()` 48532fbc7fSSatish Balay @*/ 49532fbc7fSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) 509371c9d4SSatish Balay PetscBool PetscIsInfReal(PetscReal a) { 519f4f8022SLisandro Dalcin return isinfq(a) ? PETSC_TRUE : PETSC_FALSE; 52532fbc7fSSatish Balay } 539f4f8022SLisandro Dalcin #elif defined(PETSC_HAVE_ISINF) 549371c9d4SSatish Balay PetscBool PetscIsInfReal(PetscReal a) { 559f4f8022SLisandro Dalcin return isinf(a) ? PETSC_TRUE : PETSC_FALSE; 56532fbc7fSSatish Balay } 579f4f8022SLisandro Dalcin #elif defined(PETSC_HAVE__FINITE) 58532fbc7fSSatish Balay #if defined(PETSC_HAVE_FLOAT_H) 59aaa7dc30SBarry Smith #include <float.h> /* Microsoft Windows defines _finite() in float.h */ 60532fbc7fSSatish Balay #endif 61532fbc7fSSatish Balay #if defined(PETSC_HAVE_IEEEFP_H) 62aaa7dc30SBarry Smith #include <ieeefp.h> /* Solaris prototypes these here */ 63532fbc7fSSatish Balay #endif 649371c9d4SSatish Balay PetscBool PetscIsInfReal(PetscReal a) { 659f4f8022SLisandro Dalcin return !_finite(a) ? PETSC_TRUE : PETSC_FALSE; 66532fbc7fSSatish Balay } 67532fbc7fSSatish Balay #else 689371c9d4SSatish Balay PetscBool PetscIsInfReal(PetscReal a) { 699f4f8022SLisandro Dalcin return (a && a / 2 == a) ? PETSC_TRUE : PETSC_FALSE; 70532fbc7fSSatish Balay } 71532fbc7fSSatish Balay #endif 72532fbc7fSSatish Balay 73bae46576SBarry Smith /*@C 743948c36eSLisandro Dalcin PetscIsNanReal - Returns whether the input is a Not-a-Number (NaN) value. 75bae46576SBarry Smith 76bae46576SBarry Smith Input Parameter: 77bae46576SBarry Smith . a - the floating point number 78bae46576SBarry Smith 79*811af0c4SBarry Smith Developer Notes: 80*811af0c4SBarry Smith Uses the C99 standard `isnan()` on systems where it exists. 81*811af0c4SBarry Smith 82*811af0c4SBarry Smith Otherwise uses (a != a), note that some optimizing compilers compile 83bae46576SBarry Smith out this form, thus removing the check. 84bae46576SBarry Smith 85bae46576SBarry Smith Level: beginner 86*811af0c4SBarry Smith 87*811af0c4SBarry Smith .seealso: `PetscIsNormalReal()`, `PetscIsInfReal()` 88bae46576SBarry Smith @*/ 89bae46576SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 909371c9d4SSatish Balay PetscBool PetscIsNanReal(PetscReal a) { 913948c36eSLisandro Dalcin return isnanq(a) ? PETSC_TRUE : PETSC_FALSE; 92bae46576SBarry Smith } 933948c36eSLisandro Dalcin #elif defined(PETSC_HAVE_ISNAN) 949371c9d4SSatish Balay PetscBool PetscIsNanReal(PetscReal a) { 953948c36eSLisandro Dalcin return isnan(a) ? PETSC_TRUE : PETSC_FALSE; 96bae46576SBarry Smith } 973948c36eSLisandro Dalcin #elif defined(PETSC_HAVE__ISNAN) 98bae46576SBarry Smith #if defined(PETSC_HAVE_FLOAT_H) 993948c36eSLisandro Dalcin #include <float.h> /* Microsoft Windows defines _isnan() in float.h */ 100bae46576SBarry Smith #endif 101bae46576SBarry Smith #if defined(PETSC_HAVE_IEEEFP_H) 102bae46576SBarry Smith #include <ieeefp.h> /* Solaris prototypes these here */ 103bae46576SBarry Smith #endif 1049371c9d4SSatish Balay PetscBool PetscIsNanReal(PetscReal a) { 1053948c36eSLisandro Dalcin return _isnan(a) ? PETSC_TRUE : PETSC_FALSE; 106bae46576SBarry Smith } 107bae46576SBarry Smith #else 1089371c9d4SSatish Balay PetscBool PetscIsNanReal(PetscReal a) { 1093948c36eSLisandro Dalcin return (a != a) ? PETSC_TRUE : PETSC_FALSE; 110bae46576SBarry Smith } 111bae46576SBarry Smith #endif 112