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 #if KOKKOS_VERSION >= 40799
15 #include "KokkosKernels_ArithTraits.hpp"
16 #else
17 #include "Kokkos_ArithTraits.hpp"
18 #endif
19 #include "KokkosBatched_Vector.hpp"
20 
21 //----------------------------------------------------------------------------
22 // Specializations of Kokkos{Kernels}::ArithTraits for Sacado::UQ::PCE scalar type
23 //----------------------------------------------------------------------------
24 
25 #if KOKKOS_VERSION >= 40799
26 namespace KokkosKernels {
27 #else
28 namespace Kokkos {
29 #endif
30 
31 template <typename S>
32 class ArithTraits< Sacado::UQ::PCE<S> > {
33 public:
35 
38  typedef ArithTraits<base_value_type> BAT;
39 
40  typedef typename BAT::mag_type mag_type;
41  //typedef val_type mag_type;
42 
43  static const bool is_specialized = true;
44  static const bool is_signed = BAT::is_signed;
45  static const bool is_integer = BAT::is_integer;
46  static const bool is_exact = BAT::is_exact;
47  static const bool is_complex = BAT::is_complex;
48 
49  static KOKKOS_FORCEINLINE_FUNCTION bool isInf (const val_type& x) {
50  bool res = false;
51  for (ordinal_type i=0; i<x.size(); ++i)
52  res = res || BAT::isInf(x.fastAccessCoeff(i));
53  return res;
54  }
55  static KOKKOS_FORCEINLINE_FUNCTION bool isNan (const val_type& x) {
56  bool res = false;
57  for (ordinal_type i=0; i<x.size(); ++i)
58  res = res || BAT::isInf(x.fastAccessCoeff(i));
59  return res;
60  }
61  static KOKKOS_FORCEINLINE_FUNCTION mag_type abs (const val_type& x) {
62  //return std::abs(x);
63  const ordinal_type sz = x.size();
64  mag_type n = mag_type(0);
65  for (ordinal_type i=0; i<sz; ++i)
66  n += BAT::abs( x.fastAccessCoeff(i) );
67  return n;
68  }
69  static KOKKOS_FORCEINLINE_FUNCTION val_type zero () {
70  return val_type(0.0);
71  }
72  static KOKKOS_FORCEINLINE_FUNCTION val_type one () {
73  return val_type(1.0);
74  }
75  static KOKKOS_FORCEINLINE_FUNCTION val_type min () {
76  return BAT::min();
77  }
78  static KOKKOS_FORCEINLINE_FUNCTION val_type max () {
79  return BAT::max();
80  }
81  static KOKKOS_FORCEINLINE_FUNCTION val_type real (const val_type& x) {
82  const ordinal_type sz = x.size();
83  val_type y(x.cijk(), sz);
84  for (ordinal_type i=0; i<sz; ++i)
85  y.fastAccessCoeff(i) = BAT::real(x.fastAccessCoeff(i));
86  return y;
87  }
88  static KOKKOS_FORCEINLINE_FUNCTION val_type imag (const val_type& x) {
89  const ordinal_type sz = x.size();
90  val_type y(x.cijk(), sz);
91  for (ordinal_type i=0; i<sz; ++i)
92  y.fastAccessCoeff(i) = BAT::imag(x.fastAccessCoeff(i));
93  return y;
94  }
95  static KOKKOS_FORCEINLINE_FUNCTION val_type conj (const val_type& x) {
96  const ordinal_type sz = x.size();
97  val_type y(x.cijk(), sz);
98  for (ordinal_type i=0; i<sz; ++i)
99  y.fastAccessCoeff(i) = BAT::conj(x.fastAccessCoeff(i));
100  return y;
101  }
102  static KOKKOS_FORCEINLINE_FUNCTION val_type pow (const val_type& x,
103  const val_type& y) {
104  return std::pow(x, y);
105  }
106  static KOKKOS_FORCEINLINE_FUNCTION val_type sqrt (const val_type& x) {
107  return std::sqrt(x);
108  }
109  static KOKKOS_FORCEINLINE_FUNCTION val_type log (const val_type& x) {
110  return std::log(x);
111  }
112  static KOKKOS_FORCEINLINE_FUNCTION val_type log10 (const val_type& x) {
113  return std::log10(x);
114  }
115  static KOKKOS_FORCEINLINE_FUNCTION val_type nan () {
116  return BAT::nan();
117  }
118  static KOKKOS_FORCEINLINE_FUNCTION mag_type epsilon () {
119  return BAT::epsilon();
120  }
121 
122  // Backwards compatibility with Teuchos::ScalarTraits.
124  typedef typename BAT::halfPrecision base_half_precision;
125  typedef typename BAT::doublePrecision base_double_precision;
126  typedef typename Sacado::mpl::apply<S,ordinal_type,base_half_precision>::type half_storage;
127  typedef typename Sacado::mpl::apply<S,ordinal_type,base_double_precision>::type double_storage;
130  static const bool isComplex = is_complex;
131  static const bool isOrdinal = is_integer;
132  static const bool isComparable = BAT::isComparable;
133  static const bool hasMachineParameters = BAT::hasMachineParameters;
134  static bool isnaninf (const val_type& x) {
135  return isNan (x) || isInf (x);
136  }
137  static KOKKOS_FORCEINLINE_FUNCTION mag_type magnitude (const val_type& x) {
138  return abs (x);
139  }
140  static KOKKOS_FORCEINLINE_FUNCTION val_type conjugate (const val_type& x) {
141  return conj (x);
142  }
143  static std::string name () {
144  return Sacado::StringName<val_type>::eval();
145  }
146  static KOKKOS_FORCEINLINE_FUNCTION val_type squareroot (const val_type& x) {
147  return sqrt (x);
148  }
149  static KOKKOS_FORCEINLINE_FUNCTION mag_type eps () {
150  return epsilon ();
151  }
152  static KOKKOS_FORCEINLINE_FUNCTION mag_type sfmin () {
153  return BAT::sfmin();
154  }
155  static KOKKOS_FORCEINLINE_FUNCTION int base () {
156  return BAT::base();
157  }
158  static KOKKOS_FORCEINLINE_FUNCTION mag_type prec () {
159  return BAT::prec();
160  }
161  static KOKKOS_FORCEINLINE_FUNCTION int t () {
162  return BAT::t();
163  }
164  static KOKKOS_FORCEINLINE_FUNCTION mag_type rnd () {
165  return BAT::rnd();
166  }
167  static KOKKOS_FORCEINLINE_FUNCTION int emin () {
168  return BAT::emin();
169  }
170  static KOKKOS_FORCEINLINE_FUNCTION mag_type rmin () {
171  return BAT::rmin();
172  }
173  static KOKKOS_FORCEINLINE_FUNCTION int emax () {
174  return BAT::emax();
175  }
176  static KOKKOS_FORCEINLINE_FUNCTION mag_type rmax () {
177  return BAT::rmax();
178  }
179 };
180 
181 } // namespace Kokkos
182 
183 namespace KokkosBatched {
184 
185  template <typename S>
186  struct MagnitudeScalarType< Sacado::UQ::PCE<S> > {
188 #if KOKKOS_VERSION >= 40799
189  typedef typename KokkosKernels::ArithTraits<val_type>::mag_type type;
190 #else
191  typedef typename Kokkos::ArithTraits<val_type>::mag_type type;
192 #endif
193  };
194 
195 }
196 
197 #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()