49 #ifndef _TEUCHOS_SCALARTRAITS_HPP_
50 #define _TEUCHOS_SCALARTRAITS_HPP_
58 #ifdef HAVE_TEUCHOSCORE_KOKKOSCORE
59 #include "Kokkos_Complex.hpp"
60 #endif // HAVE_TEUCHOSCORE_KOKKOSCORE
62 #ifdef HAVE_TEUCHOS_ARPREC
63 #include <arprec/mp_real.h>
66 #ifdef HAVE_TEUCHOSCORE_QUADMATH
84 operator<< (std::ostream& out,
const __float128& x);
100 #endif // HAVE_TEUCHOSCORE_QUADMATH
102 #ifdef HAVE_TEUCHOS_QD
103 #include <qd/qd_real.h>
104 #include <qd/dd_real.h>
107 #ifdef HAVE_TEUCHOS_GNU_MP
119 #ifndef DOXYGEN_SHOULD_SKIP_THIS
123 void throwScalarTraitsNanInfError(
const std::string &errMsg );
126 template<
class Scalar>
127 bool generic_real_isnaninf(
const Scalar &x)
129 #ifdef HAVE_TEUCHOSCORE_CXX11
130 if (std::isnan(x))
return true;
131 if (std::isinf(x))
return true;
134 typedef std::numeric_limits<Scalar> STD_NL;
136 const Scalar tol = 1.0;
137 if (!(x <= tol) && !(x > tol))
return true;
139 Scalar z =
static_cast<Scalar
>(0.0) * x;
140 if (!(z <= tol) && !(z > tol))
return true;
142 if (x == STD_NL::infinity() || x == -STD_NL::infinity())
return true;
149 #define TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR( VALUE, MSG ) \
150 if (isnaninf(VALUE)) { \
151 std::ostringstream omsg; \
153 Teuchos::throwScalarTraitsNanInfError(omsg.str()); \
158 struct ScalarTraits<char>
169 static inline magnitudeType magnitude(
char a) {
return static_cast<char>(std::fabs(static_cast<double>(a))); }
170 static inline char zero() {
return 0; }
171 static inline char one() {
return 1; }
172 static inline char conjugate(
char x) {
return x; }
173 static inline char real(
char x) {
return x; }
174 static inline char imag(
char) {
return 0; }
175 static inline bool isnaninf(
char ) {
return false; }
176 static inline void seedrandom(
unsigned int s) {
185 static inline char random() {
return std::rand(); }
186 static inline std::string
name() {
return "char"; }
187 static inline char squareroot(
char x) {
return (
char) std::sqrt((
double) x); }
188 static inline char pow(
char x,
char y) {
return (
char) std::pow((
double)x,(
double)y); }
189 static inline char log(
char x) {
return static_cast<char> (std::log (static_cast<double> (x))); }
190 static inline char log10(
char x) {
return static_cast<char> (std::log10 (static_cast<double> (x))); }
195 struct ScalarTraits<short int>
206 static inline magnitudeType magnitude(
short int a) {
return static_cast<short int>(std::fabs(static_cast<double>(a))); }
207 static inline short int zero() {
return 0; }
208 static inline short int one() {
return 1; }
209 static inline short int conjugate(
short int x) {
return x; }
210 static inline short int real(
short int x) {
return x; }
211 static inline short int imag(
short int) {
return 0; }
212 static inline bool isnaninf(
short int) {
return false; }
213 static inline void seedrandom(
unsigned int s) {
222 static inline short int random() {
return std::rand(); }
223 static inline std::string
name() {
return "short int"; }
224 static inline short int squareroot(
short int x) {
return (
short int) std::sqrt((
double) x); }
225 static inline short int pow(
short int x,
short int y) {
return (
short int) std::pow((
double)x,(
double)y); }
226 static inline short int log(
short int x) {
return static_cast<short int> (std::log (static_cast<double> (x))); }
227 static inline short int log10(
short int x) {
return static_cast<short int> (std::log10 (static_cast<double> (x))); }
231 struct ScalarTraits<unsigned short int>
242 static inline magnitudeType magnitude(
unsigned short int a) {
return static_cast<unsigned short int>(std::fabs(static_cast<double>(a))); }
243 static inline unsigned short int zero() {
return 0; }
244 static inline unsigned short int one() {
return 1; }
245 static inline unsigned short int conjugate(
unsigned short int x) {
return x; }
246 static inline unsigned short int real(
unsigned short int x) {
return x; }
247 static inline unsigned short int imag(
unsigned short int) {
return 0; }
248 static inline bool isnaninf(
unsigned short int) {
return false; }
249 static inline void seedrandom(
unsigned int s) {
258 static inline unsigned short int random() {
return std::rand(); }
259 static inline std::string
name() {
return "unsigned short int"; }
260 static inline unsigned short int squareroot(
unsigned short int x) {
return (
unsigned short int) std::sqrt((
double) x); }
261 static inline unsigned short int pow(
unsigned short int x,
unsigned short int y) {
return (
unsigned short int) std::pow((
double)x,(
double)y); }
262 static inline unsigned short int log(
unsigned short int x) {
return static_cast<unsigned short int> (std::log (static_cast<double> (x))); }
263 static inline unsigned short int log10(
unsigned short int x) {
return static_cast<unsigned short int> (std::log10 (static_cast<double> (x))); }
268 struct ScalarTraits<int>
279 static inline magnitudeType magnitude(
int a) {
return static_cast<int>(std::fabs(static_cast<double>(a))); }
280 static inline int zero() {
return 0; }
281 static inline int one() {
return 1; }
282 static inline int conjugate(
int x) {
return x; }
283 static inline int real(
int x) {
return x; }
284 static inline int imag(
int) {
return 0; }
285 static inline bool isnaninf(
int) {
return false; }
286 static inline void seedrandom(
unsigned int s) {
295 static inline int random() {
return std::rand(); }
296 static inline std::string
name() {
return "int"; }
297 static inline int squareroot(
int x) {
return (
int) std::sqrt((
double) x); }
298 static inline int pow(
int x,
int y) {
return (
int) std::pow((
double)x,(
double)y); }
299 static inline int log(
int x) {
return static_cast<int> (std::log (static_cast<double> (x))); }
300 static inline int log10(
int x) {
return static_cast<int> (std::log10 (static_cast<double> (x))); }
305 struct ScalarTraits<unsigned int>
316 static inline magnitudeType magnitude(
unsigned int a) {
return static_cast<unsigned int>(std::fabs(static_cast<double>(a))); }
317 static inline unsigned int zero() {
return 0; }
318 static inline unsigned int one() {
return 1; }
319 static inline unsigned int conjugate(
unsigned int x) {
return x; }
320 static inline unsigned int real(
unsigned int x) {
return x; }
321 static inline unsigned int imag(
unsigned int) {
return 0; }
322 static inline bool isnaninf(
unsigned int) {
return false; }
323 static inline void seedrandom(
unsigned int s) {
332 static inline unsigned int random() {
return std::rand(); }
333 static inline std::string
name() {
return "unsigned int"; }
334 static inline unsigned int squareroot(
unsigned int x) {
return (
unsigned int) std::sqrt((
double) x); }
335 static inline unsigned int pow(
unsigned int x,
unsigned int y) {
return (
unsigned int) std::pow((
double)x,(
double)y); }
336 static inline unsigned int log(
unsigned int x) {
return static_cast<unsigned int> (std::log (static_cast<double> (x))); }
337 static inline unsigned int log10(
unsigned int x) {
return static_cast<unsigned int> (std::log10 (static_cast<double> (x))); }
342 struct ScalarTraits<long int>
353 static inline magnitudeType magnitude(
long int a) {
return static_cast<long int>(std::fabs(static_cast<double>(a))); }
354 static inline long int zero() {
return 0; }
355 static inline long int one() {
return 1; }
356 static inline long int conjugate(
long int x) {
return x; }
357 static inline long int real(
long int x) {
return x; }
358 static inline long int imag(
long int) {
return 0; }
359 static inline bool isnaninf(
long int) {
return false; }
360 static inline void seedrandom(
unsigned int s) {
369 static inline long int random() {
return std::rand(); }
370 static inline std::string
name() {
return "long int"; }
371 static inline long int squareroot(
long int x) {
return (
long int) std::sqrt((
double) x); }
372 static inline long int pow(
long int x,
long int y) {
return (
long int) std::pow((
double)x,(
double)y); }
375 static inline long int log(
long int x) {
return static_cast<long int> (std::log (static_cast<double> (x))); }
376 static inline long int log10(
long int x) {
return static_cast<long int> (std::log10 (static_cast<double> (x))); }
381 struct ScalarTraits<long unsigned int>
392 static inline magnitudeType magnitude(
long unsigned int a) {
return static_cast<long unsigned int>(std::fabs(static_cast<double>(a))); }
393 static inline long unsigned int zero() {
return 0; }
394 static inline long unsigned int one() {
return 1; }
395 static inline long unsigned int conjugate(
long unsigned int x) {
return x; }
396 static inline long unsigned int real(
long unsigned int x) {
return x; }
397 static inline long unsigned int imag(
long unsigned int) {
return 0; }
398 static inline bool isnaninf(
long unsigned int) {
return false; }
399 static inline void seedrandom(
unsigned int s) {
408 static inline long unsigned int random() {
return std::rand(); }
409 static inline std::string
name() {
return "long unsigned int"; }
410 static inline long unsigned int squareroot(
long unsigned int x) {
return (
long unsigned int) std::sqrt((
double) x); }
411 static inline long unsigned int pow(
long unsigned int x,
long unsigned int y) {
return (
long unsigned int) std::pow((
double)x,(
double)y); }
414 static inline long unsigned int log(
long unsigned int x) {
return static_cast<long unsigned int> (std::log (static_cast<double> (x))); }
415 static inline long unsigned int log10(
long unsigned int x) {
return static_cast<long unsigned int> (std::log10 (static_cast<double> (x))); }
420 struct ScalarTraits<long long int>
431 static inline magnitudeType magnitude(
long long int a) {
return static_cast<long long int>(std::fabs(static_cast<double>(a))); }
432 static inline long long int zero() {
return 0; }
433 static inline long long int one() {
return 1; }
434 static inline long long int conjugate(
long long int x) {
return x; }
435 static inline long long int real(
long long int x) {
return x; }
436 static inline long long int imag(
long long int) {
return 0; }
437 static inline bool isnaninf(
long long int) {
return false; }
438 static inline void seedrandom(
unsigned int s) {
447 static inline long long int random() {
return std::rand(); }
448 static inline std::string
name() {
return "long long int"; }
449 static inline long long int squareroot(
long long int x) {
return (
long long int) std::sqrt((
double) x); }
450 static inline long long int pow(
long long int x,
long long int y) {
return (
long long int) std::pow((
double)x,(
double)y); }
453 static inline long long int log(
long long int x) {
return static_cast<long long int> (std::log (static_cast<double> (x))); }
454 static inline long long int log10(
long long int x) {
return static_cast<long long int> (std::log10 (static_cast<double> (x))); }
458 struct ScalarTraits<unsigned long long int>
469 static inline magnitudeType magnitude(
unsigned long long int a) {
return static_cast<unsigned long long int>(std::fabs(static_cast<double>(a))); }
470 static inline unsigned long long int zero() {
return 0; }
471 static inline unsigned long long int one() {
return 1; }
472 static inline unsigned long long int conjugate(
unsigned long long int x) {
return x; }
473 static inline unsigned long long int real(
unsigned long long int x) {
return x; }
474 static inline unsigned long long int imag(
unsigned long long int) {
return 0; }
475 static inline bool isnaninf(
unsigned long long int) {
return false; }
476 static inline void seedrandom(
unsigned int s) {
485 static inline unsigned long long int random() {
return std::rand(); }
486 static inline std::string
name() {
return "unsigned long long int"; }
487 static inline unsigned long long int squareroot(
unsigned long long int x) {
return (
unsigned long long int) std::sqrt((
double) x); }
488 static inline unsigned long long int pow(
unsigned long long int x,
unsigned long long int y) {
return (
unsigned long long int) std::pow((
double)x,(
double)y); }
491 static inline unsigned long long int log(
unsigned long long int x) {
return static_cast<unsigned long long int> (std::log (static_cast<double> (x))); }
492 static inline unsigned long long int log10(
unsigned long long int x) {
return static_cast<unsigned long long int> (std::log10 (static_cast<double> (x))); }
496 #ifdef HAVE_TEUCHOS___INT64
499 struct ScalarTraits<__int64>
510 static inline magnitudeType magnitude(__int64 a) {
return static_cast<__int64
>(std::fabs(static_cast<double>(a))); }
511 static inline __int64
zero() {
return 0; }
512 static inline __int64
one() {
return 1; }
513 static inline __int64
conjugate(__int64 x) {
return x; }
514 static inline __int64
real(__int64 x) {
return x; }
515 static inline __int64
imag(__int64) {
return 0; }
516 static inline void seedrandom(
unsigned int s) {
525 static inline __int64
random() {
return std::rand(); }
526 static inline std::string
name() {
return "__int64"; }
527 static inline __int64
squareroot(__int64 x) {
return (__int64) std::sqrt((
double) x); }
528 static inline __int64
pow(__int64 x, __int64 y) {
return (__int64) std::pow((
double)x,(
double)y); }
531 static inline __int64 log(__int64 x) {
return static_cast<__int64
> (std::log (static_cast<double> (x))); }
532 static inline __int64 log10(__int64 x) {
return static_cast<__int64
> (std::log10 (static_cast<double> (x))); }
536 struct ScalarTraits<unsigned __int64>
547 static inline magnitudeType magnitude(
unsigned __int64 a) {
return static_cast<unsigned __int64
>(std::fabs(static_cast<double>(a))); }
548 static inline unsigned __int64
zero() {
return 0; }
549 static inline unsigned __int64
one() {
return 1; }
550 static inline unsigned __int64
conjugate(
unsigned __int64 x) {
return x; }
551 static inline unsigned __int64
real(
unsigned __int64 x) {
return x; }
552 static inline unsigned __int64
imag(
unsigned __int64) {
return 0; }
553 static inline void seedrandom(
unsigned int s) {
562 static inline unsigned __int64
random() {
return std::rand(); }
563 static inline std::string
name() {
return "unsigned __int64"; }
564 static inline unsigned __int64
squareroot(
unsigned __int64 x) {
return (
unsigned __int64) std::sqrt((
double) x); }
565 static inline unsigned __int64
pow(
unsigned __int64 x,
unsigned __int64 y) {
return (
unsigned __int64) std::pow((
double)x,(
double)y); }
568 static inline unsigned __int64 log(
unsigned __int64 x) {
return static_cast<unsigned __int64
> (std::log (static_cast<double> (x))); }
569 static inline unsigned __int64 log10(
unsigned __int64 x) {
return static_cast<unsigned __int64
> (std::log10 (static_cast<double> (x))); }
572 #endif // HAVE_TEUCHOS___INT64
581 struct ScalarTraits<float>
591 static inline float eps() {
592 return std::numeric_limits<float>::epsilon();
594 static inline float sfmin() {
595 return std::numeric_limits<float>::min();
597 static inline float base() {
598 return static_cast<float>(std::numeric_limits<float>::radix);
600 static inline float prec() {
603 static inline float t() {
604 return static_cast<float>(std::numeric_limits<float>::digits);
606 static inline float rnd() {
607 return ( std::numeric_limits<float>::round_style == std::round_to_nearest ?
one() :
zero() );
609 static inline float emin() {
610 return static_cast<float>(std::numeric_limits<float>::min_exponent);
612 static inline float rmin() {
613 return std::numeric_limits<float>::min();
615 static inline float emax() {
616 return static_cast<float>(std::numeric_limits<float>::max_exponent);
618 static inline float rmax() {
619 return std::numeric_limits<float>::max();
624 TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR(
625 a,
"Error, the input value to magnitude(...) a = " << a <<
" can not be NaN!" );
629 static inline float zero() {
return(0.0
f); }
630 static inline float one() {
return(1.0
f); }
631 static inline float conjugate(
float x) {
return(x); }
632 static inline float real(
float x) {
return x; }
633 static inline float imag(
float) {
return zero(); }
634 static inline float nan() {
636 return 0.0f/std::sin(0.0
f);
638 return std::numeric_limits<float>::quiet_NaN();
641 static inline bool isnaninf(
float x) {
642 return generic_real_isnaninf<float>(x);
644 static inline void seedrandom(
unsigned int s) {
652 static inline float random() {
float rnd = (float) std::rand() / RAND_MAX;
return (-1.0
f + 2.0
f * rnd); }
653 static inline std::string
name() {
return "float"; }
657 TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR(
658 x,
"Error, the input value to squareroot(...) x = " << x <<
" can not be NaN!" );
661 const float rtn = std::sqrt(x);
666 static inline float pow(
float x,
float y) {
return std::pow(x,y); }
667 static inline float pi() {
return 3.14159265358979323846f; }
668 static inline float log(
float x) {
return std::log(x); }
669 static inline float log10(
float x) {
return std::log10(x); }
679 struct ScalarTraits<double>
692 #if defined(HAVE_TEUCHOS_DOUBLE_TO_QD)
694 #elif defined(HAVE_TEUCHOS_DOUBLE_TO_ARPREC)
704 static inline double eps() {
705 return std::numeric_limits<double>::epsilon();
707 static inline double sfmin() {
708 return std::numeric_limits<double>::min();
710 static inline double base() {
711 return std::numeric_limits<double>::radix;
713 static inline double prec() {
716 static inline double t() {
717 return std::numeric_limits<double>::digits;
719 static inline double rnd() {
720 return ( std::numeric_limits<double>::round_style == std::round_to_nearest ?
double(1.0) :
double(0.0) );
722 static inline double emin() {
723 return std::numeric_limits<double>::min_exponent;
725 static inline double rmin() {
726 return std::numeric_limits<double>::min();
728 static inline double emax() {
729 return std::numeric_limits<double>::max_exponent;
731 static inline double rmax() {
732 return std::numeric_limits<double>::max();
737 TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR(
738 a,
"Error, the input value to magnitude(...) a = " << a <<
" can not be NaN!" );
742 static inline double zero() {
return 0.0; }
743 static inline double one() {
return 1.0; }
744 static inline double conjugate(
double x) {
return(x); }
745 static inline double real(
double x) {
return(x); }
746 static inline double imag(
double) {
return(0); }
747 static inline double nan() {
749 return 0.0/std::sin(0.0);
751 return std::numeric_limits<double>::quiet_NaN();
754 static inline bool isnaninf(
double x) {
755 return generic_real_isnaninf<double>(x);
757 static inline void seedrandom(
unsigned int s) {
765 static inline double random() {
double rnd = (double) std::rand() / RAND_MAX;
return (
double)(-1.0 + 2.0 *
rnd); }
766 static inline std::string
name() {
return "double"; }
770 TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR(
771 x,
"Error, the input value to squareroot(...) x = " << x <<
" can not be NaN!" );
774 const double rtn = std::sqrt(x);
779 static inline double pow(
double x,
double y) {
return std::pow(x,y); }
780 static inline double pi() {
return 3.14159265358979323846; }
781 static inline double log(
double x) {
return std::log(x); }
782 static inline double log10(
double x) {
return std::log10(x); }
786 #ifdef HAVE_TEUCHOSCORE_QUADMATH
789 struct ScalarTraits<__float128> {
803 static __float128
eps () {
804 return FLT128_EPSILON;
806 static __float128
sfmin () {
809 static __float128
base () {
812 static __float128
prec () {
815 static __float128
t () {
816 return FLT128_MANT_DIG;
818 static __float128
rnd () {
821 static __float128
emin () {
822 return FLT128_MIN_EXP;
824 static __float128
rmin () {
827 static __float128
emax () {
828 return FLT128_MAX_EXP;
830 static __float128
rmax () {
836 static __float128
zero () {
839 static __float128
one () {
842 static __float128
conjugate (
const __float128& x) {
845 static __float128
real (
const __float128& x) {
848 static __float128
imag (
const __float128& ) {
851 static __float128
nan () {
852 return strtoflt128 (
"NAN()", NULL);
854 static bool isnaninf (
const __float128& x) {
855 return isinfq (x) || isnanq (x);
857 static inline void seedrandom (
unsigned int s) {
865 static __float128
random () {
868 const __float128 scalingFactor =
869 static_cast<__float128
> (std::numeric_limits<double>::min ()) /
870 static_cast<__float128> (2.0);
871 const __float128 higherOrderTerm =
873 const __float128 lowerOrderTerm =
876 return higherOrderTerm + lowerOrderTerm;
878 static std::string
name () {
881 static __float128
squareroot (
const __float128& x) {
884 static __float128
pow (
const __float128& x,
const __float128& y) {
887 static __float128
pi() {
return 3.14159265358979323846; }
888 static __float128 log (
const __float128& x) {
891 static __float128 log10 (
const __float128& x) {
895 #endif // HAVE_TEUCHOSCORE_QUADMATH
899 #ifdef HAVE_TEUCHOS_QD
901 bool operator&&(
const dd_real &a,
const dd_real &b);
902 bool operator&&(
const qd_real &a,
const qd_real &b);
905 struct ScalarTraits<dd_real>
915 static inline dd_real
eps() {
return std::numeric_limits<dd_real>::epsilon(); }
916 static inline dd_real
sfmin() {
return std::numeric_limits<dd_real>::min(); }
917 static inline dd_real
base() {
return std::numeric_limits<dd_real>::radix; }
918 static inline dd_real
prec() {
return eps()*
base(); }
919 static inline dd_real
t() {
return std::numeric_limits<dd_real>::digits; }
920 static inline dd_real
rnd() {
return ( std::numeric_limits<dd_real>::round_style == std::round_to_nearest ? dd_real(1.0) : dd_real(0.0) ); }
921 static inline dd_real
emin() {
return std::numeric_limits<dd_real>::min_exponent; }
922 static inline dd_real
rmin() {
return std::numeric_limits<dd_real>::min(); }
923 static inline dd_real
emax() {
return std::numeric_limits<dd_real>::max_exponent; }
924 static inline dd_real
rmax() {
return std::numeric_limits<dd_real>::max(); }
928 TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR(
929 a,
"Error, the input value to magnitude(...) a = " << a <<
" can not be NaN!" );
933 static inline dd_real
zero() {
return dd_real(0.0); }
934 static inline dd_real
one() {
return dd_real(1.0); }
935 static inline dd_real
conjugate(dd_real x) {
return(x); }
936 static inline dd_real
real(dd_real x) {
return x ; }
937 static inline dd_real
imag(dd_real) {
return zero(); }
938 static inline dd_real
nan() {
return dd_real::_nan; }
939 static inline bool isnaninf(dd_real x) {
return isnan(x) || isinf(x); }
940 static inline void seedrandom(
unsigned int s) {
949 static inline dd_real
random() {
return ddrand(); }
950 static inline std::string
name() {
return "dd_real"; }
954 TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR(
955 x,
"Error, the input value to squareroot(...) x = " << x <<
" can not be NaN!" );
959 static inline dd_real
pow(dd_real x, dd_real y) { return ::pow(x,y); }
960 static inline dd_real
pi() {
return 3.14159265358979323846; }
962 static inline dd_real log(dd_real x) { return ::log(x); }
963 static inline dd_real log10(dd_real x) { return ::log10(x); }
968 struct ScalarTraits<qd_real>
978 static inline qd_real
eps() {
return std::numeric_limits<qd_real>::epsilon(); }
979 static inline qd_real
sfmin() {
return std::numeric_limits<qd_real>::min(); }
980 static inline qd_real
base() {
return std::numeric_limits<qd_real>::radix; }
981 static inline qd_real
prec() {
return eps()*
base(); }
982 static inline qd_real
t() {
return std::numeric_limits<qd_real>::digits; }
983 static inline qd_real
rnd() {
return ( std::numeric_limits<qd_real>::round_style == std::round_to_nearest ? qd_real(1.0) : qd_real(0.0) ); }
984 static inline qd_real
emin() {
return std::numeric_limits<qd_real>::min_exponent; }
985 static inline qd_real
rmin() {
return std::numeric_limits<qd_real>::min(); }
986 static inline qd_real
emax() {
return std::numeric_limits<qd_real>::max_exponent; }
987 static inline qd_real
rmax() {
return std::numeric_limits<qd_real>::max(); }
991 TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR(
992 a,
"Error, the input value to magnitude(...) a = " << a <<
" can not be NaN!" );
996 static inline qd_real
zero() {
return qd_real(0.0); }
997 static inline qd_real
one() {
return qd_real(1.0); }
998 static inline qd_real
conjugate(qd_real x) {
return(x); }
999 static inline qd_real
real(qd_real x) {
return x ; }
1000 static inline qd_real
imag(qd_real) {
return zero(); }
1001 static inline qd_real
nan() {
return qd_real::_nan; }
1002 static inline bool isnaninf(qd_real x) {
return isnan(x) || isinf(x); }
1003 static inline void seedrandom(
unsigned int s) {
1012 static inline qd_real
random() {
return qdrand(); }
1013 static inline std::string
name() {
return "qd_real"; }
1016 #ifdef TEUCHOS_DEBUG
1017 TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR(
1018 x,
"Error, the input value to squareroot(...) x = " << x <<
" can not be NaN!" );
1022 static inline qd_real
pow(qd_real x, qd_real y) { return ::pow(x,y); }
1023 static inline qd_real
pi() {
return 3.14159265358979323846; }
1025 static inline qd_real log(qd_real x) { return ::log(x); }
1026 static inline qd_real log10(qd_real x) { return ::log10(x); }
1030 #endif // HAVE_TEUCHOS_QD
1033 #ifdef HAVE_TEUCHOS_GNU_MP
1036 extern gmp_randclass gmp_rng;
1056 struct ScalarTraits<mpf_class>
1066 static inline mpf_class
zero() { mpf_class
zero = 0.0;
return zero; }
1067 static inline mpf_class
one() { mpf_class
one = 1.0;
return one; }
1068 static inline mpf_class
conjugate(mpf_class x) {
return x; }
1069 static inline mpf_class
real(mpf_class x) {
return(x); }
1070 static inline mpf_class
imag(mpf_class x) {
return(0); }
1071 static inline bool isnaninf(mpf_class x) {
return false; }
1072 static inline void seedrandom(
unsigned int s) {
1073 unsigned long int seedVal =
static_cast<unsigned long int>(s);
1074 gmp_rng.seed( seedVal );
1076 static inline mpf_class
random() {
1077 return gmp_rng.get_f();
1079 static inline std::string
name() {
return "mpf_class"; }
1080 static inline mpf_class
squareroot(mpf_class x) {
return std::sqrt(x); }
1081 static inline mpf_class
pow(mpf_class x, mpf_class y) {
return pow(x,y); }
1085 #endif // HAVE_TEUCHOS_GNU_MP
1087 #ifdef HAVE_TEUCHOS_ARPREC
1092 struct ScalarTraits<mp_real>
1104 static inline mp_real
zero() { mp_real
zero = 0.0;
return zero; }
1105 static inline mp_real
one() { mp_real
one = 1.0;
return one; }
1106 static inline mp_real
conjugate(mp_real x) {
return x; }
1107 static inline mp_real
real(mp_real x) {
return(x); }
1108 static inline mp_real
imag(mp_real x) {
return zero(); }
1109 static inline bool isnaninf(mp_real x) {
return false; }
1110 static inline void seedrandom(
unsigned int s) {
1111 long int seedVal =
static_cast<long int>(s);
1114 static inline mp_real
random() {
return mp_rand(); }
1115 static inline std::string
name() {
return "mp_real"; }
1116 static inline mp_real
squareroot(mp_real x) {
return sqrt(x); }
1117 static inline mp_real
pow(mp_real x, mp_real y) {
return pow(x,y); }
1118 static inline mp_real
pi() {
return 3.14159265358979323846; }
1123 #endif // HAVE_TEUCHOS_ARPREC
1126 #ifdef HAVE_TEUCHOS_COMPLEX
1131 struct ScalarTraits<
1135 typedef std::complex<T> ComplexT;
1136 typedef std::complex<typename ScalarTraits<T>::halfPrecision>
halfPrecision;
1137 typedef std::complex<typename ScalarTraits<T>::doublePrecision>
doublePrecision;
1138 typedef typename ScalarTraits<T>::magnitudeType
magnitudeType;
1141 static const bool isOrdinal = ScalarTraits<T>::isOrdinal;
1156 #ifdef TEUCHOS_DEBUG
1157 TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR(
1158 a,
"Error, the input value to magnitude(...) a = " << a <<
" can not be NaN!" );
1162 static inline ComplexT
zero() {
return ComplexT(ScalarTraits<magnitudeType>::zero(),ScalarTraits<magnitudeType>::zero()); }
1163 static inline ComplexT
one() {
return ComplexT(ScalarTraits<magnitudeType>::one(),ScalarTraits<magnitudeType>::zero()); }
1164 static inline ComplexT
conjugate(ComplexT a){
return ComplexT(a.real(),-a.imag()); }
1167 static inline ComplexT
nan() {
return ComplexT(ScalarTraits<magnitudeType>::nan(),ScalarTraits<magnitudeType>::nan()); }
1170 static inline ComplexT
random()
1174 return ComplexT(rnd1,rnd2);
1176 static inline std::string
name() {
return std::string(
"std::complex<")+std::string(ScalarTraits<magnitudeType>::name())+std::string(
">"); }
1178 static inline ComplexT
squareroot(ComplexT x)
1180 #ifdef TEUCHOS_DEBUG
1181 TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR(
1182 x,
"Error, the input value to squareroot(...) x = " << x <<
" can not be NaN!" );
1184 typedef ScalarTraits<magnitudeType> STMT;
1185 const T r = x.real(), i = x.imag(),
zero = STMT::zero(), two = 2.0;
1186 const T a = STMT::squareroot((r*r)+(i*i));
1187 const T nr = STMT::squareroot((a+r)/two);
1188 const T ni = ( i ==
zero ?
zero : STMT::squareroot((a-r)/two) );
1189 return ComplexT(nr,ni);
1203 static inline ComplexT
pow(ComplexT x, ComplexT y) {
return pow(x,y); }
1206 #endif // HAVE_TEUCHOS_COMPLEX
1208 #ifdef HAVE_TEUCHOSCORE_KOKKOSCORE
1211 struct ScalarTraits<
1215 typedef Kokkos::complex<T> ComplexT;
1216 typedef Kokkos::complex<typename ScalarTraits<T>::halfPrecision>
halfPrecision;
1217 typedef Kokkos::complex<typename ScalarTraits<T>::doublePrecision>
doublePrecision;
1218 typedef typename ScalarTraits<T>::magnitudeType
magnitudeType;
1221 static const bool isOrdinal = ScalarTraits<T>::isOrdinal;
1236 #ifdef TEUCHOS_DEBUG
1237 TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR(
1238 a,
"Error, the input value to magnitude(...) a = " << a <<
" can not be NaN!" );
1240 return std::abs(std::complex<T>(a));
1242 static inline ComplexT
zero() {
return ComplexT(ScalarTraits<magnitudeType>::zero(),ScalarTraits<magnitudeType>::zero()); }
1243 static inline ComplexT
one() {
return ComplexT(ScalarTraits<magnitudeType>::one(),ScalarTraits<magnitudeType>::zero()); }
1244 static inline ComplexT
conjugate(ComplexT a){
return ComplexT(a.real(),-a.imag()); }
1247 static inline ComplexT
nan() {
return ComplexT(ScalarTraits<magnitudeType>::nan(),ScalarTraits<magnitudeType>::nan()); }
1250 static inline ComplexT
random()
1254 return ComplexT(rnd1,rnd2);
1256 static inline std::string
name() {
return std::string(
"Kokkos::complex<")+std::string(ScalarTraits<magnitudeType>::name())+std::string(
">"); }
1258 static inline ComplexT
squareroot(ComplexT x)
1260 #ifdef TEUCHOS_DEBUG
1261 TEUCHOS_SCALAR_TRAITS_NAN_INF_ERR(
1262 x,
"Error, the input value to squareroot(...) x = " << x <<
" can not be NaN!" );
1264 typedef ScalarTraits<magnitudeType> STMT;
1265 const T r = x.real(), i = x.imag(),
zero = STMT::zero(), two = 2.0;
1266 const T a = STMT::squareroot((r*r)+(i*i));
1267 const T nr = STMT::squareroot((a+r)/two);
1268 const T ni = ( i ==
zero ?
zero : STMT::squareroot((a-r)/two) );
1269 return ComplexT(nr,ni);
1271 static inline ComplexT
pow(ComplexT x, ComplexT y) {
return pow(std::complex<T>(x), std::complex<T>(y)); }
1274 #endif // HAVE_TEUCHOSCORE_KOKKOSCORE
1276 #endif // DOXYGEN_SHOULD_SKIP_THIS
1280 #endif // _TEUCHOS_SCALARTRAITS_HPP_
T magnitudeType
Mandatory typedef for result of magnitude.
static magnitudeType eps()
Returns relative machine precision.
static T squareroot(T x)
Returns a number of magnitudeType that is the square root of this scalar type x.
static const bool isComparable
Determines if scalar type supports relational operators such as <, >, <=, >=.
static magnitudeType real(T a)
Returns the real part of the scalar type a.
Teuchos header file which uses auto-configuration information to include necessary C++ headers...
Declaration and default implementation for basic traits for the scalar field type.
static T pow(T x, T y)
Returns the result of raising one scalar x to the power y.
static magnitudeType emax()
Returns the largest exponent before overflow.
static magnitudeType base()
Returns the base of the machine.
static const bool hasMachineParameters
Determines if scalar type have machine-specific parameters (i.e. eps(), sfmin(), base(), prec(), t(), rnd(), emin(), rmin(), emax(), rmax() are supported).
ostream & operator<<(ostream &os, const pair< Packet, Packet > &arg)
static std::string name()
Returns the name of this scalar type.
static magnitudeType rmax()
Overflow theshold - (base^emax)*(1-eps)
static magnitudeType sfmin()
Returns safe minimum (sfmin), such that 1/sfmin does not overflow.
static T pi()
Returns the value of PI.
T coordinateType
Typedef for coordinates.
static T conjugate(T a)
Returns the conjugate of the scalar type a.
static const bool isOrdinal
Determines if scalar type is an ordinal type.
static magnitudeType prec()
Returns eps*base.
#define TEUCHOSCORE_LIB_DLL_EXPORT
static magnitudeType t()
Returns the number of (base) digits in the mantissa.
static magnitudeType rmin()
Returns the underflow threshold - base^(emin-1)
T doublePrecision
Typedef for double precision.
static void seedrandom(unsigned int s)
Seed the random number generator returned by random().
static magnitudeType imag(T a)
Returns the imaginary part of the scalar type a.
static bool isnaninf(const T &x)
Returns true if x is NaN or Inf.
static magnitudeType magnitude(T a)
Returns the magnitudeType of the scalar type a.
static T nan()
Returns a number that represents NaN.
std::istream & operator>>(std::istream &in, CustomDataType &object)
static T zero()
Returns representation of zero for this scalar type.
static T random()
Returns a random number (between -one() and +one()) of this scalar type.
T halfPrecision
Typedef for half precision.
static const bool isComplex
Determines if scalar type is std::complex.
static magnitudeType emin()
Returns the minimum exponent before (gradual) underflow.
static magnitudeType rnd()
Returns 1.0 when rounding occurs in addition, 0.0 otherwise.
static T one()
Returns representation of one for this scalar type.