Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Kokkos_ArithTraits_UQ_PCE.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef KOKKOS_ARITHTRAITS_UQ_PCE_HPP
11 #define KOKKOS_ARITHTRAITS_UQ_PCE_HPP
12 
13 #include "Sacado_UQ_PCE.hpp"
14 #include "Kokkos_ArithTraits.hpp"
15 #include "KokkosBatched_Vector.hpp"
16 
17 //----------------------------------------------------------------------------
18 // Specializations of Kokkos::ArithTraits for Sacado::UQ::PCE scalar type
19 //----------------------------------------------------------------------------
20 
21 namespace Kokkos {
22 
23 template <typename S>
24 class ArithTraits< Sacado::UQ::PCE<S> > {
25 public:
27 
30  typedef ArithTraits<base_value_type> BAT;
31 
32  typedef typename BAT::mag_type mag_type;
33  //typedef val_type mag_type;
34 
35  static const bool is_specialized = true;
36  static const bool is_signed = BAT::is_signed;
37  static const bool is_integer = BAT::is_integer;
38  static const bool is_exact = BAT::is_exact;
39  static const bool is_complex = BAT::is_complex;
40 
41  static KOKKOS_FORCEINLINE_FUNCTION bool isInf (const val_type& x) {
42  bool res = false;
43  for (ordinal_type i=0; i<x.size(); ++i)
44  res = res || BAT::isInf(x.fastAccessCoeff(i));
45  return res;
46  }
47  static KOKKOS_FORCEINLINE_FUNCTION bool isNan (const val_type& x) {
48  bool res = false;
49  for (ordinal_type i=0; i<x.size(); ++i)
50  res = res || BAT::isInf(x.fastAccessCoeff(i));
51  return res;
52  }
53  static KOKKOS_FORCEINLINE_FUNCTION mag_type abs (const val_type& x) {
54  //return std::abs(x);
55  const ordinal_type sz = x.size();
56  mag_type n = mag_type(0);
57  for (ordinal_type i=0; i<sz; ++i)
58  n += BAT::abs( x.fastAccessCoeff(i) );
59  return n;
60  }
61  static KOKKOS_FORCEINLINE_FUNCTION val_type zero () {
62  return val_type(0.0);
63  }
64  static KOKKOS_FORCEINLINE_FUNCTION val_type one () {
65  return val_type(1.0);
66  }
67  static KOKKOS_FORCEINLINE_FUNCTION val_type min () {
68  return BAT::min();
69  }
70  static KOKKOS_FORCEINLINE_FUNCTION val_type max () {
71  return BAT::max();
72  }
73  static KOKKOS_FORCEINLINE_FUNCTION val_type real (const val_type& x) {
74  const ordinal_type sz = x.size();
75  val_type y(x.cijk(), sz);
76  for (ordinal_type i=0; i<sz; ++i)
77  y.fastAccessCoeff(i) = BAT::real(x.fastAccessCoeff(i));
78  return y;
79  }
80  static KOKKOS_FORCEINLINE_FUNCTION val_type imag (const val_type& x) {
81  const ordinal_type sz = x.size();
82  val_type y(x.cijk(), sz);
83  for (ordinal_type i=0; i<sz; ++i)
84  y.fastAccessCoeff(i) = BAT::imag(x.fastAccessCoeff(i));
85  return y;
86  }
87  static KOKKOS_FORCEINLINE_FUNCTION val_type conj (const val_type& x) {
88  const ordinal_type sz = x.size();
89  val_type y(x.cijk(), sz);
90  for (ordinal_type i=0; i<sz; ++i)
91  y.fastAccessCoeff(i) = BAT::conj(x.fastAccessCoeff(i));
92  return y;
93  }
94  static KOKKOS_FORCEINLINE_FUNCTION val_type pow (const val_type& x,
95  const val_type& y) {
96  return std::pow(x, y);
97  }
98  static KOKKOS_FORCEINLINE_FUNCTION val_type sqrt (const val_type& x) {
99  return std::sqrt(x);
100  }
101  static KOKKOS_FORCEINLINE_FUNCTION val_type log (const val_type& x) {
102  return std::log(x);
103  }
104  static KOKKOS_FORCEINLINE_FUNCTION val_type log10 (const val_type& x) {
105  return std::log10(x);
106  }
107  static KOKKOS_FORCEINLINE_FUNCTION val_type nan () {
108  return BAT::nan();
109  }
110  static KOKKOS_FORCEINLINE_FUNCTION mag_type epsilon () {
111  return BAT::epsilon();
112  }
113 
114  // Backwards compatibility with Teuchos::ScalarTraits.
116  typedef typename BAT::halfPrecision base_half_precision;
117  typedef typename BAT::doublePrecision base_double_precision;
118  typedef typename Sacado::mpl::apply<S,ordinal_type,base_half_precision>::type half_storage;
119  typedef typename Sacado::mpl::apply<S,ordinal_type,base_double_precision>::type double_storage;
122  static const bool isComplex = is_complex;
123  static const bool isOrdinal = is_integer;
124  static const bool isComparable = BAT::isComparable;
125  static const bool hasMachineParameters = BAT::hasMachineParameters;
126  static bool isnaninf (const val_type& x) {
127  return isNan (x) || isInf (x);
128  }
129  static KOKKOS_FORCEINLINE_FUNCTION mag_type magnitude (const val_type& x) {
130  return abs (x);
131  }
132  static KOKKOS_FORCEINLINE_FUNCTION val_type conjugate (const val_type& x) {
133  return conj (x);
134  }
135  static std::string name () {
136  return Sacado::StringName<val_type>::eval();
137  }
138  static KOKKOS_FORCEINLINE_FUNCTION val_type squareroot (const val_type& x) {
139  return sqrt (x);
140  }
141  static KOKKOS_FORCEINLINE_FUNCTION mag_type eps () {
142  return epsilon ();
143  }
144  static KOKKOS_FORCEINLINE_FUNCTION mag_type sfmin () {
145  return BAT::sfmin();
146  }
147  static KOKKOS_FORCEINLINE_FUNCTION int base () {
148  return BAT::base();
149  }
150  static KOKKOS_FORCEINLINE_FUNCTION mag_type prec () {
151  return BAT::prec();
152  }
153  static KOKKOS_FORCEINLINE_FUNCTION int t () {
154  return BAT::t();
155  }
156  static KOKKOS_FORCEINLINE_FUNCTION mag_type rnd () {
157  return BAT::rnd();
158  }
159  static KOKKOS_FORCEINLINE_FUNCTION int emin () {
160  return BAT::emin();
161  }
162  static KOKKOS_FORCEINLINE_FUNCTION mag_type rmin () {
163  return BAT::rmin();
164  }
165  static KOKKOS_FORCEINLINE_FUNCTION int emax () {
166  return BAT::emax();
167  }
168  static KOKKOS_FORCEINLINE_FUNCTION mag_type rmax () {
169  return BAT::rmax();
170  }
171 };
172 
173 } // namespace Kokkos
174 
175 namespace KokkosBatched {
176 
177  template <typename S>
178  struct MagnitudeScalarType< Sacado::UQ::PCE<S> > {
180  typedef typename Kokkos::ArithTraits<val_type>::mag_type type;
181  };
182 
183 }
184 
185 #endif /* #ifndef KOKKOS_ARITHTRAITS_UQ_PCE_HPP */
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
static KOKKOS_FORCEINLINE_FUNCTION val_type max()
static KOKKOS_FORCEINLINE_FUNCTION val_type min()
static KOKKOS_FORCEINLINE_FUNCTION mag_type rmax()
static KOKKOS_FORCEINLINE_FUNCTION mag_type abs(const val_type &x)
static KOKKOS_FORCEINLINE_FUNCTION val_type log10(const val_type &x)
KOKKOS_INLINE_FUNCTION PCE< Storage > pow(const PCE< Storage > &a, const PCE< Storage > &b)
static KOKKOS_FORCEINLINE_FUNCTION mag_type rmin()
static KOKKOS_FORCEINLINE_FUNCTION val_type one()
static KOKKOS_FORCEINLINE_FUNCTION int base()
Sacado::mpl::apply< S, ordinal_type, base_double_precision >::type double_storage
static KOKKOS_FORCEINLINE_FUNCTION val_type log(const val_type &x)
static KOKKOS_FORCEINLINE_FUNCTION mag_type magnitude(const val_type &x)
static KOKKOS_FORCEINLINE_FUNCTION int emax()
KOKKOS_INLINE_FUNCTION PCE< Storage > min(const typename PCE< Storage >::value_type &a, const PCE< Storage > &b)
static KOKKOS_FORCEINLINE_FUNCTION val_type squareroot(const val_type &x)
static KOKKOS_FORCEINLINE_FUNCTION int t()
static KOKKOS_FORCEINLINE_FUNCTION val_type nan()
static KOKKOS_FORCEINLINE_FUNCTION int emin()
static KOKKOS_FORCEINLINE_FUNCTION mag_type prec()
static KOKKOS_FORCEINLINE_FUNCTION val_type imag(const val_type &x)
static KOKKOS_FORCEINLINE_FUNCTION mag_type sfmin()
KOKKOS_INLINE_FUNCTION PCE< Storage > max(const typename PCE< Storage >::value_type &a, const PCE< Storage > &b)
KOKKOS_INLINE_FUNCTION PCE< Storage > abs(const PCE< Storage > &a)
static KOKKOS_FORCEINLINE_FUNCTION mag_type rnd()
static KOKKOS_FORCEINLINE_FUNCTION bool isNan(const val_type &x)
static KOKKOS_FORCEINLINE_FUNCTION val_type real(const val_type &x)
Sacado::Random< double > rnd
Sacado::mpl::apply< S, ordinal_type, base_half_precision >::type half_storage
static KOKKOS_FORCEINLINE_FUNCTION val_type conjugate(const val_type &x)
static KOKKOS_FORCEINLINE_FUNCTION val_type pow(const val_type &x, const val_type &y)
static KOKKOS_FORCEINLINE_FUNCTION mag_type eps()
static KOKKOS_FORCEINLINE_FUNCTION bool isInf(const val_type &x)
static KOKKOS_FORCEINLINE_FUNCTION val_type conj(const val_type &x)
KOKKOS_INLINE_FUNCTION PCE< Storage > log(const PCE< Storage > &a)
static KOKKOS_FORCEINLINE_FUNCTION val_type zero()
KOKKOS_INLINE_FUNCTION PCE< Storage > log10(const PCE< Storage > &a)
int n
static KOKKOS_FORCEINLINE_FUNCTION val_type sqrt(const val_type &x)
static KOKKOS_FORCEINLINE_FUNCTION mag_type epsilon()