42 #ifndef KOKKOS_ARITHTRAITS_UQ_PCE_HPP 
   43 #define KOKKOS_ARITHTRAITS_UQ_PCE_HPP 
   46 #include "Kokkos_ArithTraits.hpp" 
   47 #include "KokkosBatched_Vector.hpp" 
   57 class ArithTraits< Sacado::UQ::PCE<S> > {
 
   63   typedef ArithTraits<base_value_type> 
BAT;
 
   68   static const bool is_specialized = 
true;
 
   69   static const bool is_signed = BAT::is_signed;
 
   70   static const bool is_integer = BAT::is_integer;
 
   71   static const bool is_exact = BAT::is_exact;
 
   72   static const bool is_complex = BAT::is_complex;
 
   77       res = res || BAT::isInf(x.fastAccessCoeff(i));
 
   83       res = res || BAT::isInf(x.fastAccessCoeff(i));
 
   91       n += 
BAT::abs( x.fastAccessCoeff(i) );
 
  110       y.fastAccessCoeff(i) = BAT::real(x.fastAccessCoeff(i));
 
  117       y.fastAccessCoeff(i) = BAT::imag(x.fastAccessCoeff(i));
 
  124       y.fastAccessCoeff(i) = BAT::conj(x.fastAccessCoeff(i));
 
  144     return BAT::epsilon();
 
  151   typedef typename Sacado::mpl::apply<S,ordinal_type,base_half_precision>::type 
half_storage;
 
  152   typedef typename Sacado::mpl::apply<S,ordinal_type,base_double_precision>::type 
double_storage;
 
  155   static const bool isComplex = is_complex;
 
  156   static const bool isOrdinal = is_integer;
 
  157   static const bool isComparable = BAT::isComparable;
 
  158   static const bool hasMachineParameters = BAT::hasMachineParameters;
 
  160     return isNan (x) || isInf (x);
 
  169     return Sacado::StringName<val_type>::eval();
 
  180   static KOKKOS_FORCEINLINE_FUNCTION 
int base () {
 
  186   static KOKKOS_FORCEINLINE_FUNCTION 
int t () {
 
  192   static KOKKOS_FORCEINLINE_FUNCTION 
int emin () {
 
  198   static KOKKOS_FORCEINLINE_FUNCTION 
int emax () {
 
  209 namespace KokkosBatched {
 
  211   template <
typename S>
 
  212   struct MagnitudeScalarType< Sacado::UQ::PCE<S> > {
 
  214     typedef typename Kokkos::Details::ArithTraits<val_type>::mag_type 
type;
 
static KOKKOS_FORCEINLINE_FUNCTION mag_type epsilon()
 
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
 
val_type::ordinal_type ordinal_type
 
static KOKKOS_FORCEINLINE_FUNCTION int emin()
 
static KOKKOS_FORCEINLINE_FUNCTION val_type log(const val_type &x)
 
Sacado::mpl::apply< S, ordinal_type, base_double_precision >::type double_storage
 
static KOKKOS_FORCEINLINE_FUNCTION int base()
 
static KOKKOS_FORCEINLINE_FUNCTION val_type max()
 
KOKKOS_INLINE_FUNCTION PCE< Storage > pow(const PCE< Storage > &a, const PCE< Storage > &b)
 
static KOKKOS_FORCEINLINE_FUNCTION mag_type magnitude(const val_type &x)
 
static KOKKOS_FORCEINLINE_FUNCTION mag_type eps()
 
static KOKKOS_FORCEINLINE_FUNCTION val_type conjugate(const val_type &x)
 
Kokkos::Details::ArithTraits< val_type >::mag_type type
 
static KOKKOS_FORCEINLINE_FUNCTION val_type sqrt(const val_type &x)
 
BAT::doublePrecision base_double_precision
 
static KOKKOS_FORCEINLINE_FUNCTION val_type min()
 
ArithTraits< base_value_type > BAT
 
KOKKOS_INLINE_FUNCTION PCE< Storage > min(const typename PCE< Storage >::value_type &a, const PCE< Storage > &b)
 
Sacado::UQ::PCE< half_storage > halfPrecision
 
static KOKKOS_FORCEINLINE_FUNCTION val_type squareroot(const val_type &x)
 
static KOKKOS_FORCEINLINE_FUNCTION val_type imag(const val_type &x)
 
KOKKOS_INLINE_FUNCTION PCE< Storage > max(const typename PCE< Storage >::value_type &a, const PCE< Storage > &b)
 
Sacado::UQ::PCE< S > val_type
 
static KOKKOS_FORCEINLINE_FUNCTION int emax()
 
static KOKKOS_FORCEINLINE_FUNCTION mag_type prec()
 
KOKKOS_INLINE_FUNCTION PCE< Storage > abs(const PCE< Storage > &a)
 
static KOKKOS_FORCEINLINE_FUNCTION val_type log10(const val_type &x)
 
static KOKKOS_FORCEINLINE_FUNCTION int t()
 
Sacado::Random< double > rnd
 
Sacado::UQ::PCE< S > val_type
 
Sacado::mpl::apply< S, ordinal_type, base_half_precision >::type half_storage
 
static KOKKOS_FORCEINLINE_FUNCTION mag_type rmin()
 
static bool isnaninf(const val_type &x)
 
static KOKKOS_FORCEINLINE_FUNCTION val_type zero()
 
static KOKKOS_FORCEINLINE_FUNCTION bool isNan(const val_type &x)
 
static KOKKOS_FORCEINLINE_FUNCTION bool isInf(const val_type &x)
 
static KOKKOS_FORCEINLINE_FUNCTION val_type pow(const val_type &x, const val_type &y)
 
BAT::halfPrecision base_half_precision
 
static KOKKOS_FORCEINLINE_FUNCTION val_type nan()
 
static std::string name()
 
static KOKKOS_FORCEINLINE_FUNCTION val_type one()
 
static KOKKOS_FORCEINLINE_FUNCTION mag_type rmax()
 
static KOKKOS_FORCEINLINE_FUNCTION mag_type sfmin()
 
static KOKKOS_FORCEINLINE_FUNCTION val_type conj(const val_type &x)
 
static KOKKOS_FORCEINLINE_FUNCTION mag_type abs(const val_type &x)
 
val_type::value_type base_value_type
 
KOKKOS_INLINE_FUNCTION PCE< Storage > log(const PCE< Storage > &a)
 
static KOKKOS_FORCEINLINE_FUNCTION mag_type rnd()
 
KOKKOS_INLINE_FUNCTION PCE< Storage > log10(const PCE< Storage > &a)
 
static KOKKOS_FORCEINLINE_FUNCTION val_type real(const val_type &x)
 
Sacado::UQ::PCE< double_storage > doublePrecision