xref: /petsc/src/sys/utils/mathinf.c (revision 8b49ba18373ad29c29af9cb0a122b33db25f4ce7)
1532fbc7fSSatish Balay #include <petscsys.h>
2532fbc7fSSatish Balay /*@C
3*8b49ba18SBarry Smith       PetscIsNormal - Returns PETSC_TRUE if the input value satisfies isnormal()
4532fbc7fSSatish Balay 
5532fbc7fSSatish Balay     Input Parameter:
6*8b49ba18SBarry Smith .     a - the PetscRealValue
7532fbc7fSSatish Balay 
8*8b49ba18SBarry Smith      Notes: uses the C99 standard isnormal() on systems where they exist.
9*8b49ba18SBarry Smith       Uses isnormalq() with __float128
10*8b49ba18SBarry Smith       Otherwises always returns true
11*8b49ba18SBarry Smith 
12*8b49ba18SBarry Smith      Level: beginner
13*8b49ba18SBarry Smith @*/
14*8b49ba18SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
15*8b49ba18SBarry Smith PetscBool PetscIsNormalScalar(PetscScalar a)
16*8b49ba18SBarry Smith {
17*8b49ba18SBarry Smith   return isnormalq(PetscAbsScalar(a));
18*8b49ba18SBarry Smith }
19*8b49ba18SBarry Smith PetscBool PetscIsNormalReal(PetscReal a)
20*8b49ba18SBarry Smith {
21*8b49ba18SBarry Smith   return isnormalq(a);
22*8b49ba18SBarry Smith }
23*8b49ba18SBarry Smith #elif defined(PETSC_HAVE_ISNORMAL)
24*8b49ba18SBarry Smith PetscBool PetscIsNormalScalar(PetscScalar a)
25*8b49ba18SBarry Smith {
26*8b49ba18SBarry Smith   return isnormal(PetscAbsScalar(a));
27*8b49ba18SBarry Smith }
28*8b49ba18SBarry Smith PetscBool PetscIsNormalReal(PetscReal a)
29*8b49ba18SBarry Smith {
30*8b49ba18SBarry Smith   return isnormal(a);
31*8b49ba18SBarry Smith }
32*8b49ba18SBarry Smith #else
33*8b49ba18SBarry Smith PetscBool PetscIsNormalScalar(PetscScalar a)
34*8b49ba18SBarry Smith {
35*8b49ba18SBarry Smith   return PETSC_TRUE;
36*8b49ba18SBarry Smith }
37*8b49ba18SBarry Smith PetscBool PetscIsNormalReal(PetscReal a)
38*8b49ba18SBarry Smith {
39*8b49ba18SBarry Smith   return PETSC_TRUE;
40*8b49ba18SBarry Smith }
41*8b49ba18SBarry Smith #endif
42*8b49ba18SBarry Smith 
43*8b49ba18SBarry Smith /*@C
44*8b49ba18SBarry Smith       PetscIsInfOrNan - Returns an error code if the input double has an infinity for Not-a-number (Nan) value, otherwise 0.
45*8b49ba18SBarry Smith 
46*8b49ba18SBarry Smith     Input Parameter:
47*8b49ba18SBarry Smith .     a - the floating point number
48532fbc7fSSatish Balay 
49532fbc7fSSatish Balay      Notes: uses the C99 standard isinf() and isnan() on systems where they exist.
50532fbc7fSSatish Balay       Otherwises uses ((a - a) != 0.0), note that some optimizing compiles compile
51532fbc7fSSatish Balay       out this form, thus removing the check.
52532fbc7fSSatish Balay 
53532fbc7fSSatish Balay      Level: beginner
54532fbc7fSSatish Balay @*/
55532fbc7fSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128)
56f48dab2eSKarl Rupp PetscErrorCode PetscIsInfOrNanScalar(PetscScalar a)
57f48dab2eSKarl Rupp {
58532fbc7fSSatish Balay   return isinfq(PetscAbsScalar(a)) || isnanq(PetscAbsScalar(a));
59532fbc7fSSatish Balay }
60f48dab2eSKarl Rupp PetscErrorCode PetscIsInfOrNanReal(PetscReal a)
61f48dab2eSKarl Rupp {
62532fbc7fSSatish Balay   return isinfq(a) || isnanq(a);
63532fbc7fSSatish Balay }
64532fbc7fSSatish Balay #elif defined(PETSC_HAVE_ISINF) && defined(PETSC_HAVE_ISNAN)
65f48dab2eSKarl Rupp PetscErrorCode PetscIsInfOrNanScalar(PetscScalar a)
66f48dab2eSKarl Rupp {
67532fbc7fSSatish Balay   return isinf(PetscAbsScalar(a)) || isnan(PetscAbsScalar(a));
68532fbc7fSSatish Balay }
69f48dab2eSKarl Rupp PetscErrorCode PetscIsInfOrNanReal(PetscReal a)
70f48dab2eSKarl Rupp {
71532fbc7fSSatish Balay   return isinf(a) || isnan(a);
72532fbc7fSSatish Balay }
73532fbc7fSSatish Balay #elif defined(PETSC_HAVE__FINITE) && defined(PETSC_HAVE__ISNAN)
74532fbc7fSSatish Balay #if defined(PETSC_HAVE_FLOAT_H)
75532fbc7fSSatish Balay #include "float.h"  /* Microsoft Windows defines _finite() in float.h */
76532fbc7fSSatish Balay #endif
77532fbc7fSSatish Balay #if defined(PETSC_HAVE_IEEEFP_H)
78532fbc7fSSatish Balay #include "ieeefp.h"  /* Solaris prototypes these here */
79532fbc7fSSatish Balay #endif
80f48dab2eSKarl Rupp PetscErrorCode PetscIsInfOrNanScalar(PetscScalar a)
81f48dab2eSKarl Rupp {
82532fbc7fSSatish Balay   return !_finite(PetscAbsScalar(a)) || _isnan(PetscAbsScalar(a));
83532fbc7fSSatish Balay }
84f48dab2eSKarl Rupp PetscErrorCode PetscIsInfOrNanReal(PetscReal a)
85f48dab2eSKarl Rupp {
86532fbc7fSSatish Balay   return !_finite(a) || _isnan(a);
87532fbc7fSSatish Balay }
88532fbc7fSSatish Balay #else
89f48dab2eSKarl Rupp PetscErrorCode PetscIsInfOrNanScalar(PetscScalar a)
90f48dab2eSKarl Rupp {
91532fbc7fSSatish Balay   return  ((a - a) != (PetscScalar)0);
92532fbc7fSSatish Balay }
93f48dab2eSKarl Rupp PetscErrorCode PetscIsInfOrNanReal(PetscReal a)
94f48dab2eSKarl Rupp {
95532fbc7fSSatish Balay   return ((a - a) != 0);
96532fbc7fSSatish Balay }
97532fbc7fSSatish Balay #endif
98532fbc7fSSatish Balay 
99