42 #ifndef SACADO_MP_SCALAR_TRAITS_IMP_HPP
43 #define SACADO_MP_SCALAR_TRAITS_IMP_HPP
50 #include "Sacado_mpl_apply.hpp"
55 template <
typename S,
bool reduct_across_vector>
72 typedef typename Sacado::mpl::apply<S,ordinal_type,value_mag_type>::type
storage_mag_type;
73 typedef typename Sacado::mpl::apply<S,ordinal_type,value_half_type>::type
storage_half_type;
81 static const bool isComplex = TVT::isComplex;
82 static const bool isOrdinal = TVT::isOrdinal;
83 static const bool isComparable = TVT::isComparable;
84 static const bool hasMachineParameters = TVT::hasMachineParameters;
124 for (
int i=0; i<sz; i++)
125 y.fastAccessCoeff(i) = TVT::conjugate(x.fastAccessCoeff(i));
154 for (
int i=0; i<x.size(); i++)
155 if (TVT::isnaninf(x.fastAccessCoeff(i)))
160 static void seedrandom(
unsigned int s) { TVT::seedrandom(s); }
164 static const char *
name() {
return "Sacado::MP::Vector<>"; }
181 template <
typename S>
201 static const bool isComplex = TVT::isComplex;
202 static const bool isOrdinal = TVT::isOrdinal;
203 static const bool isComparable = TVT::isComparable;
204 static const bool hasMachineParameters = TVT::hasMachineParameters;
238 for (
int i=0; i<sz; i++)
239 y.fastAccessCoeff(i) = TVT::conjugate(x.fastAccessCoeff(i));
246 for (
int i=0; i<sz; i++)
247 y.fastAccessCoeff(i) = TVT::real(x.fastAccessCoeff(i));
254 for (
int i=0; i<sz; i++)
255 y.fastAccessCoeff(i) = TVT::imag(x.fastAccessCoeff(i));
262 for (
int i=0; i<x.size(); i++)
263 if (TVT::isnaninf(x.fastAccessCoeff(i)))
268 static void seedrandom(
unsigned int s) { TVT::seedrandom(s); }
272 static const char *
name() {
return "Sacado::MP::Vector<>"; }
287 template <
typename Ordinal,
typename VecType,
typename Serializer>
293 typedef typename Sacado::ValueType<VecType>::type
ValueT;
312 const VecType buffer[],
317 for (
Ordinal i=0; i<count; i++) {
318 int my_sz = buffer[i].size();
320 if (sz == 0) tot_sz = my_sz;
322 if (tot_sz != my_sz) {
331 Ordinal b2 = vs.fromCountToIndirectBytes(tot_sz, cx->coeff());
343 const VecType buffer[],
349 for (
Ordinal i=0; i<count; i++) {
351 int my_sz = buffer[i].size();
353 if (sz == 0) tot_sz = my_sz;
359 if (tot_sz != my_sz) {
368 Ordinal b2 = vs.fromCountToIndirectBytes(tot_sz, cx->coeff());
372 vs.serialize(tot_sz, cx->coeff(), b2, charBuffer);
382 const char charBuffer[],
386 while (bytes_used < bytes) {
409 const char charBuffer[],
413 for (
Ordinal i=0; i<count; i++) {
422 if (sz == 0) tot_sz = *my_sz;
423 buffer[i] = VecType(tot_sz,
ValueT(0.0));
429 vs.deserialize(*b2, charBuffer, *my_sz, buffer[i].coeff());
439 template <
typename Ordinal,
typename VecType,
bool is_static = false>
445 typedef typename Sacado::ValueType<VecType>::type
ValueT;
467 const VecType buffer[]) {
474 const VecType buffer[],
483 const char charBuffer[]) {
490 const char charBuffer[],
502 template <
typename Ordinal,
typename VecType>
504 typedef typename Sacado::ValueType<VecType>::type
ValueT;
511 vSerT::supportsDirectSerialization;
518 return DSerT::fromCountToDirectBytes(count);
523 return DSerT::convertToCharPtr(ptr);
528 return DSerT::convertToCharPtr(ptr);
533 return DSerT::fromDirectBytesToCount(bytes);
538 return DSerT::convertFromCharPtr(ptr);
543 return DSerT::convertFromCharPtr(ptr);
553 const VecType buffer[]) {
555 return DSerT::fromCountToIndirectBytes(count, buffer);
557 return STI::fromCountToIndirectBytes(count, buffer);
562 const VecType buffer[],
566 return DSerT::serialize(count, buffer, bytes, charBuffer);
568 return STI::serialize(count, buffer, bytes, charBuffer);
573 const char charBuffer[]) {
575 return DSerT::fromIndirectBytesToCount(bytes, charBuffer);
577 return STI::fromIndirectBytesToCount(bytes, charBuffer);
582 const char charBuffer[],
586 return DSerT::deserialize(bytes, charBuffer, count, buffer);
588 return STI::deserialize(bytes, charBuffer, count, buffer);
596 template <
typename Ordinal,
typename VecType,
typename ValueSerializer>
636 const VecType buffer[])
const {
642 const VecType buffer[],
644 char charBuffer[])
const {
650 const char charBuffer[])
const {
656 const char charBuffer[],
658 VecType buffer[])
const {
670 #endif // SACADO_MP_SCALAR_TRAITS_IMP_HPP
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > fabs(const PCE< Storage > &a)
SerializationImp< Ordinal, VecType, ValueSerializer > Imp
Implementation.
static magnitudeType imag(const ScalarType &x)
static Ordinal fromIndirectBytesToCount(const Serializer &vs, const Ordinal bytes, const char charBuffer[], const Ordinal sz=0)
Return the number of objects for bytes of storage.
static void serialize(const Ordinal count, const intbuffer[], const Ordinal bytes, char charBuffer[])
static const VecType * convertFromCharPtr(const char *ptr)
Convert the pointer type from char*.
static value_mag_type rmax()
TVT::magnitudeType value_mag_type
Teuchos::RCP< const ValueSerializer > vs
Serializer for value types.
static char * convertToCharPtr(VecType *ptr)
Convert the pointer type to char*.
static void deserialize(const Ordinal bytes, const char charBuffer[], const Ordinal count, VecType buffer[])
Deserialize from an indirect char[] buffer.
static value_mag_type emax()
static const char * name()
static value_mag_type sfmin()
static Ordinal fromCountToIndirectBytes(const Ordinal count, const intbuffer[])
static Ordinal fromCountToDirectBytes(const Ordinal count)
static ScalarType log10(const ScalarType &x)
static void deserialize(const Ordinal bytes, const char charBuffer[], const Ordinal count, VecType buffer[])
Deserialize from an indirect char[] buffer.
static value_mag_type rmin()
static void deserialize(const Serializer &vs, const Ordinal bytes, const char charBuffer[], const Ordinal count, VecType buffer[], const Ordinal sz=0)
Deserialize from an indirect char[] buffer.
Teuchos::DirectSerializationTraits< Ordinal, VecType > DSerT
static ScalarType log10(const ScalarType &x)
TVT::doublePrecision value_double_type
static Ordinal fromCountToDirectBytes(const Ordinal count)
Return the number of bytes for count objects.
static value_mag_type rmin()
static ScalarType log(const ScalarType &x)
Teuchos::ScalarTraits< value_type >::coordinateType coordinateType
static magnitudeType magnitude(const ScalarType &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > pow(const PCE< Storage > &a, const PCE< Storage > &b)
static Ordinal fromCountToIndirectBytes(const Ordinal count, const T buffer[])
value_mag_type magnitudeType
SerializationImp< Ordinal, VecType, ValueSerializer > Imp
Implementation.
static VecType * convertFromCharPtr(char *ptr)
Convert the pointer type from char*.
static value_mag_type sfmin()
TVT::halfPrecision value_half_type
Sacado::ValueType< VecType >::type ValueT
Value type.
static value_mag_type eps()
static bool isnaninf(const ScalarType &x)
static ScalarType squareroot(const ScalarType &x)
static value_mag_type emin()
static value_mag_type prec()
static value_mag_type rnd()
Teuchos::SerializationTraits< Ordinal, int > iSerT
How to serialize ints.
static void seedrandom(unsigned int s)
Ordinal fromCountToIndirectBytes(const Ordinal count, const VecType buffer[]) const
Return the number of bytes for count objects.
static const char * name()
TVT::magnitudeType value_mag_type
Teuchos::SerializationTraits< Ordinal, ValueT > vSerT
static value_mag_type rnd()
DS::DefaultSerializerType ValueSerializer
Default serializer type for values.
static value_mag_type base()
Sacado::mpl::apply< S, ordinal_type, value_mag_type >::type storage_mag_type
static ScalarType squareroot(const ScalarType &x)
static Ordinal fromCountToIndirectBytes(const Ordinal count, const VecType buffer[])
Return the number of bytes for count objects.
static magnitudeType real(const ScalarType &x)
static void seedrandom(unsigned int s)
static DefaultSerializerType getDefaultSerializer()
static void serialize(const Serializer &vs, const Ordinal count, const VecType buffer[], const Ordinal bytes, char charBuffer[], const Ordinal sz=0)
Serialize to an indirect char[] buffer.
static value_mag_type emin()
Sacado::MP::Vector< storage_double_type > doublePrecision
Ordinal sz
Specified number of derivative components;.
SerializerImp(const Teuchos::RCP< const ValueSerializer > &vs_, Ordinal sz_=0)
Constructor.
static int * convertFromCharPtr(char *ptr)
static ScalarType conjugate(const ScalarType &x)
static ScalarType random()
static value_mag_type base()
Sacado::MP::Vector< storage_double_type > doublePrecision
Serialization implementation for all Vector types.
Ordinal getSerializerSize() const
Return specified serializer size.
Sacado::MP::Vector< S > ScalarType
static Ordinal fromDirectBytesToCount(const Ordinal bytes)
Return the number of objects for bytes of storage.
static ScalarType pow(const ScalarType &x, const ScalarType &y)
static T * convertFromCharPtr(char *ptr)
Ordinal fromIndirectBytesToCount(const Ordinal bytes, const char charBuffer[]) const
Return the number of objects for bytes of storage.
static Ordinal fromIndirectBytesToCount(const Ordinal bytes, const char charBuffer[])
Return the number of objects for bytes of storage.
TVT::doublePrecision value_double_type
Sacado::MP::Vector< storage_half_type > halfPrecision
static ScalarType random()
static Ordinal fromCountToIndirectBytes(const Ordinal count, const VecType buffer[])
Return the number of bytes for count objects.
value_type doublePrecision
static ScalarType conjugate(const ScalarType &x)
static const bool supportsDirectSerialization
Whether the type T supports direct serialization.
Teuchos::SerializationTraits< Ordinal, Ordinal > oSerT
How to serialize ordinals.
static const bool supportsDirectSerialization
Whether the type T supports direct serialization.
Sacado::MP::Vector< S > ScalarType
static ScalarType real(const ScalarType &x)
Sacado::Random< double > rnd
void serialize(const Ordinal count, const VecType buffer[], const Ordinal bytes, char charBuffer[]) const
Serialize to an indirect char[] buffer.
S::ordinal_type ordinal_type
Teuchos::ScalarTraits< value_type > TVT
static ScalarType pow(const ScalarType &x, const ScalarType &y)
Sacado::ValueType< VecType >::type ValueT
Value type of Vec type.
Implementation of Teuchos::SerializationTraits for all Vector types.
static magnitudeType magnitude(const ScalarType &a)
static const char * convertToCharPtr(const VecType *ptr)
Convert the pointer type to const char*.
static Ordinal fromCountToIndirectBytes(const Serializer &vs, const Ordinal count, const VecType buffer[], const Ordinal sz=0)
Return the number of bytes for count objects.
static Ordinal fromIndirectBytesToCount(const Ordinal bytes, const char charBuffer[])
Return the number of objects for bytes of storage.
Sacado::mpl::apply< S, ordinal_type, value_half_type >::type storage_half_type
static ScalarType log(const ScalarType &x)
Sacado::mpl::apply< S, ordinal_type, value_half_type >::type storage_half_type
TVT::halfPrecision value_half_type
Sacado::mpl::apply< S, ordinal_type, value_double_type >::type storage_double_type
static value_mag_type rmax()
static ScalarType imag(const ScalarType &x)
Teuchos::ScalarTraits< value_type >::coordinateType coordinateType
Sacado::ValueType< VecType >::type ValueT
void deserialize(const Ordinal bytes, const char charBuffer[], const Ordinal count, VecType buffer[]) const
Deserialize from an indirect char[] buffer.
static value_mag_type t()
Teuchos::DefaultSerializer< Ordinal, ValueT > DS
Default serializer for values.
An indirect serialization object for all Vector types.
Sacado::MP::Vector< storage_half_type > halfPrecision
Sacado::MP::Vector< storage_mag_type > magnitudeType
ValueSerializer value_serializer_type
Typename of value serializer.
Teuchos::RCP< const value_serializer_type > getValueSerializer() const
Get nested value serializer.
KOKKOS_INLINE_FUNCTION PCE< Storage > log(const PCE< Storage > &a)
static value_mag_type eps()
static value_mag_type t()
static const bool supportsDirectSerialization
Whether we support direct serialization.
KOKKOS_INLINE_FUNCTION PCE< Storage > log10(const PCE< Storage > &a)
Teuchos::ScalarTraits< value_type > TVT
static void serialize(const Ordinal count, const VecType buffer[], const Ordinal bytes, char charBuffer[])
Serialize to an indirect char[] buffer.
Sacado::mpl::apply< S, ordinal_type, value_double_type >::type storage_double_type
static void serialize(const Ordinal count, const T buffer[], const Ordinal bytes, char charBuffer[])
Sacado::MP::SerializationTraitsImp< Ordinal, VecType > STI
static Ordinal fromCountToDirectBytes(const Ordinal count)
static value_mag_type emax()
static bool isnaninf(const ScalarType &x)
static value_mag_type prec()
S::ordinal_type ordinal_type
Sacado::mpl::apply< S, ordinal_type, value_mag_type >::type storage_mag_type
static void serialize(const Ordinal count, const VecType buffer[], const Ordinal bytes, char charBuffer[])
Serialize to an indirect char[] buffer.