Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Sacado_MP_ScalarTraitsImp.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Stokhos Package
5 // Copyright (2009) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Eric T. Phipps (etphipp@sandia.gov).
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef SACADO_MP_SCALAR_TRAITS_IMP_HPP
43 #define SACADO_MP_SCALAR_TRAITS_IMP_HPP
44 
45 #include "Teuchos_ScalarTraits.hpp"
48 #include "Teuchos_Assert.hpp"
50 #include "Sacado_mpl_apply.hpp"
51 
52 namespace Sacado {
53  namespace MP {
54 
55  template <typename S, bool reduct_across_vector>
56  struct ScalarTraitsImp {};
57 
58  // Implementation of Teuchos::ScalarTraits where reductions are taken
59  // across the components of MP::Vector. In this case magnitudeType is
60  // a scalar
61  template <typename S>
62  struct ScalarTraitsImp<S,true> {
64  typedef typename S::value_type value_type;
65  typedef typename S::ordinal_type ordinal_type;
67 
71 
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;
74  typedef typename Sacado::mpl::apply<S,ordinal_type,value_double_type>::type storage_double_type;
75 
80 
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;
85 
86  static value_mag_type eps() { return TVT::eps(); }
87 
88  static value_mag_type sfmin() { return TVT::sfmin(); }
89 
90  static value_mag_type base() { return TVT::base(); }
91 
92  static value_mag_type prec() { return TVT::prec(); }
93 
94  static value_mag_type t() { return TVT::t(); }
95 
96  static value_mag_type rnd() { return TVT::rnd(); }
97 
98  static value_mag_type emin() { return TVT::emin(); }
99 
100  static value_mag_type rmin() { return TVT::rmin(); }
101 
102  static value_mag_type emax() { return TVT::emax(); }
103 
104  static value_mag_type rmax() { return TVT::rmax(); }
105 
106  static magnitudeType magnitude(const ScalarType& a) {
107  magnitudeType m = magnitudeType(0.0);
108  const ordinal_type sz = a.size();
109  for (ordinal_type i=0; i<sz; ++i) {
110  value_mag_type t = TVT::magnitude(a.fastAccessCoeff(i));
111  m +=t*t;
112  }
113  return std::sqrt(m);
114  }
115 
116  static ScalarType zero() { return ScalarType(0.0); }
117 
118  static ScalarType one() { return ScalarType(1.0); }
119 
120 
121  static ScalarType conjugate(const ScalarType& x) {
122  int sz = x.size();
123  ScalarType y(sz, value_type(0.0));
124  for (int i=0; i<sz; i++)
125  y.fastAccessCoeff(i) = TVT::conjugate(x.fastAccessCoeff(i));
126  return y;
127  }
128 
129 
130  static magnitudeType real(const ScalarType& x) {
131  magnitudeType m = magnitudeType(0.0);
132  const ordinal_type sz = x.size();
133  for (ordinal_type i=0; i<sz; ++i) {
134  value_mag_type t = TVT::real(x.fastAccessCoeff(i));
135  m +=t*t;
136  }
137  return std::sqrt(m);
138  }
139 
140 
141  static magnitudeType imag(const ScalarType& x) {
142  magnitudeType m = magnitudeType(0.0);
143  const ordinal_type sz = x.size();
144  for (ordinal_type i=0; i<sz; ++i) {
145  value_mag_type t = TVT::imag(x.fastAccessCoeff(i));
146  m +=t*t;
147  }
148  return std::sqrt(m);
149  }
150 
151  static value_type nan() { return TVT::nan(); }
152 
153  static bool isnaninf(const ScalarType& x) {
154  for (int i=0; i<x.size(); i++)
155  if (TVT::isnaninf(x.fastAccessCoeff(i)))
156  return true;
157  return false;
158  }
159 
160  static void seedrandom(unsigned int s) { TVT::seedrandom(s); }
161 
162  static ScalarType random() { return ScalarType(TVT::random()); }
163 
164  static const char * name() { return "Sacado::MP::Vector<>"; }
165 
166  static ScalarType squareroot(const ScalarType& x) { return std::sqrt(x); }
167 
168  static ScalarType pow(const ScalarType& x, const ScalarType& y) {
169  return std::pow(x,y);
170  }
171 
172  static ScalarType log(const ScalarType& x) { return std::log(x); }
173 
174  static ScalarType log10(const ScalarType& x) { return std::log10(x); }
175 
176  }; // class ScalarTraitsImp<S,true>
177 
178  // Implementation of Teuchos::ScalarTraits where reductions are not taken
179  // across the components of MP::Vector. In this case magnitudeType is
180  // an MP::Vector
181  template <typename S>
182  struct ScalarTraitsImp<S,false> {
184  typedef typename S::value_type value_type;
185  typedef typename S::ordinal_type ordinal_type;
187 
191 
192  typedef typename Sacado::mpl::apply<S,ordinal_type,value_mag_type>::type storage_mag_type;
193  typedef typename Sacado::mpl::apply<S,ordinal_type,value_half_type>::type storage_half_type;
194  typedef typename Sacado::mpl::apply<S,ordinal_type,value_double_type>::type storage_double_type;
195 
200 
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;
205 
206  static value_mag_type eps() { return TVT::eps(); }
207 
208  static value_mag_type sfmin() { return TVT::sfmin(); }
209 
210  static value_mag_type base() { return TVT::base(); }
211 
212  static value_mag_type prec() { return TVT::prec(); }
213 
214  static value_mag_type t() { return TVT::t(); }
215 
216  static value_mag_type rnd() { return TVT::rnd(); }
217 
218  static value_mag_type emin() { return TVT::emin(); }
219 
220  static value_mag_type rmin() { return TVT::rmin(); }
221 
222  static value_mag_type emax() { return TVT::emax(); }
223 
224  static value_mag_type rmax() { return TVT::rmax(); }
225 
226  static magnitudeType magnitude(const ScalarType& a) {
227  return std::fabs(a);
228  }
229 
230  static ScalarType zero() { return ScalarType(0.0); }
231 
232  static ScalarType one() { return ScalarType(1.0); }
233 
234 
235  static ScalarType conjugate(const ScalarType& x) {
236  int sz = x.size();
237  ScalarType y(sz, value_type(0.0));
238  for (int i=0; i<sz; i++)
239  y.fastAccessCoeff(i) = TVT::conjugate(x.fastAccessCoeff(i));
240  return y;
241  }
242 
243  static ScalarType real(const ScalarType& x) {
244  int sz = x.size();
245  ScalarType y(sz, value_type(0.0));
246  for (int i=0; i<sz; i++)
247  y.fastAccessCoeff(i) = TVT::real(x.fastAccessCoeff(i));
248  return y;
249  }
250 
251  static ScalarType imag(const ScalarType& x) {
252  int sz = x.size();
253  ScalarType y(sz, value_type(0.0));
254  for (int i=0; i<sz; i++)
255  y.fastAccessCoeff(i) = TVT::imag(x.fastAccessCoeff(i));
256  return y;
257  }
258 
259  static value_type nan() { return TVT::nan(); }
260 
261  static bool isnaninf(const ScalarType& x) {
262  for (int i=0; i<x.size(); i++)
263  if (TVT::isnaninf(x.fastAccessCoeff(i)))
264  return true;
265  return false;
266  }
267 
268  static void seedrandom(unsigned int s) { TVT::seedrandom(s); }
269 
270  static ScalarType random() { return ScalarType(TVT::random()); }
271 
272  static const char * name() { return "Sacado::MP::Vector<>"; }
273 
274  static ScalarType squareroot(const ScalarType& x) { return std::sqrt(x); }
275 
276  static ScalarType pow(const ScalarType& x, const ScalarType& y) {
277  return std::pow(x,y);
278  }
279 
280  static ScalarType log(const ScalarType& x) { return std::log(x); }
281 
282  static ScalarType log10(const ScalarType& x) { return std::log10(x); }
283 
284  }; // class ScalarTraitsImp<S,false>
285 
287  template <typename Ordinal, typename VecType, typename Serializer>
289 
290  private:
291 
293  typedef typename Sacado::ValueType<VecType>::type ValueT;
294 
297 
300 
301  public:
302 
304  static const bool supportsDirectSerialization = false;
305 
307 
308 
310  static Ordinal fromCountToIndirectBytes(const Serializer& vs,
311  const Ordinal count,
312  const VecType buffer[],
313  const Ordinal sz = 0) {
314  Ordinal bytes = 0;
315  VecType *x = NULL;
316  const VecType *cx;
317  for (Ordinal i=0; i<count; i++) {
318  int my_sz = buffer[i].size();
319  int tot_sz = sz;
320  if (sz == 0) tot_sz = my_sz;
321  Ordinal b1 = iSerT::fromCountToIndirectBytes(1, &tot_sz);
322  if (tot_sz != my_sz) {
323  if (x == NULL)
324  x = new VecType;
325  *x = buffer[i];
326  x->reset(tot_sz);
327  cx = x;
328  }
329  else
330  cx = &(buffer[i]);
331  Ordinal b2 = vs.fromCountToIndirectBytes(tot_sz, cx->coeff());
333  bytes += b1+b2+b3;
334  }
335  if (x != NULL)
336  delete x;
337  return bytes;
338  }
339 
341  static void serialize (const Serializer& vs,
342  const Ordinal count,
343  const VecType buffer[],
344  const Ordinal bytes,
345  char charBuffer[],
346  const Ordinal sz = 0) {
347  VecType *x = NULL;
348  const VecType *cx;
349  for (Ordinal i=0; i<count; i++) {
350  // First serialize size
351  int my_sz = buffer[i].size();
352  int tot_sz = sz;
353  if (sz == 0) tot_sz = my_sz;
354  Ordinal b1 = iSerT::fromCountToIndirectBytes(1, &tot_sz);
355  iSerT::serialize(1, &tot_sz, b1, charBuffer);
356  charBuffer += b1;
357 
358  // Next serialize vector coefficients
359  if (tot_sz != my_sz) {
360  if (x == NULL)
361  x = new VecType;
362  *x = buffer[i];
363  x->reset(tot_sz);
364  cx = x;
365  }
366  else
367  cx = &(buffer[i]);
368  Ordinal b2 = vs.fromCountToIndirectBytes(tot_sz, cx->coeff());
370  oSerT::serialize(1, &b2, b3, charBuffer);
371  charBuffer += b3;
372  vs.serialize(tot_sz, cx->coeff(), b2, charBuffer);
373  charBuffer += b2;
374  }
375  if (x != NULL)
376  delete x;
377  }
378 
380  static Ordinal fromIndirectBytesToCount(const Serializer& vs,
381  const Ordinal bytes,
382  const char charBuffer[],
383  const Ordinal sz = 0) {
384  Ordinal count = 0;
385  Ordinal bytes_used = 0;
386  while (bytes_used < bytes) {
387 
388  // Bytes for size
390  bytes_used += b1;
391  charBuffer += b1;
392 
393  // Bytes for vector coefficients
395  const Ordinal *b2 = oSerT::convertFromCharPtr(charBuffer);
396  bytes_used += b3;
397  charBuffer += b3;
398  bytes_used += *b2;
399  charBuffer += *b2;
400 
401  ++count;
402  }
403  return count;
404  }
405 
407  static void deserialize (const Serializer& vs,
408  const Ordinal bytes,
409  const char charBuffer[],
410  const Ordinal count,
411  VecType buffer[],
412  const Ordinal sz = 0) {
413  for (Ordinal i=0; i<count; i++) {
414 
415  // Deserialize size
417  const int *my_sz = iSerT::convertFromCharPtr(charBuffer);
418  charBuffer += b1;
419 
420  // Create empty Vector object of given size
421  int tot_sz = sz;
422  if (sz == 0) tot_sz = *my_sz;
423  buffer[i] = VecType(tot_sz, ValueT(0.0));
424 
425  // Deserialize vector coefficients
427  const Ordinal *b2 = oSerT::convertFromCharPtr(charBuffer);
428  charBuffer += b3;
429  vs.deserialize(*b2, charBuffer, *my_sz, buffer[i].coeff());
430  charBuffer += *b2;
431  }
432 
433  }
435 
436  };
437 
439  template <typename Ordinal, typename VecType, bool is_static = false>
441 
442  private:
443 
445  typedef typename Sacado::ValueType<VecType>::type ValueT;
446 
449 
452 
455 
456  public:
457 
459  static const bool supportsDirectSerialization =
461 
463 
464 
467  const VecType buffer[]) {
469  DS::getDefaultSerializer(), count, buffer);
470  }
471 
473  static void serialize (const Ordinal count,
474  const VecType buffer[],
475  const Ordinal bytes,
476  char charBuffer[]) {
478  DS::getDefaultSerializer(), count, buffer, bytes, charBuffer);
479  }
480 
483  const char charBuffer[]) {
485  DS::getDefaultSerializer(), bytes, charBuffer);
486  }
487 
489  static void deserialize (const Ordinal bytes,
490  const char charBuffer[],
491  const Ordinal count,
492  VecType buffer[]) {
494  DS::getDefaultSerializer(), bytes, charBuffer, count, buffer);
495  }
496 
498 
499  };
500 
502  template <typename Ordinal, typename VecType>
503  struct SerializationTraitsImp<Ordinal, VecType, true> {
504  typedef typename Sacado::ValueType<VecType>::type ValueT;
508 
510  static const bool supportsDirectSerialization =
511  vSerT::supportsDirectSerialization;
512 
514 
515 
517  static Ordinal fromCountToDirectBytes(const Ordinal count) {
518  return DSerT::fromCountToDirectBytes(count);
519  }
520 
522  static char* convertToCharPtr( VecType* ptr ) {
523  return DSerT::convertToCharPtr(ptr);
524  }
525 
527  static const char* convertToCharPtr( const VecType* ptr ) {
528  return DSerT::convertToCharPtr(ptr);
529  }
530 
532  static Ordinal fromDirectBytesToCount(const Ordinal bytes) {
533  return DSerT::fromDirectBytesToCount(bytes);
534  }
535 
537  static VecType* convertFromCharPtr( char* ptr ) {
538  return DSerT::convertFromCharPtr(ptr);
539  }
540 
542  static const VecType* convertFromCharPtr( const char* ptr ) {
543  return DSerT::convertFromCharPtr(ptr);
544  }
545 
547 
549 
550 
553  const VecType buffer[]) {
555  return DSerT::fromCountToIndirectBytes(count, buffer);
556  else
557  return STI::fromCountToIndirectBytes(count, buffer);
558  }
559 
561  static void serialize (const Ordinal count,
562  const VecType buffer[],
563  const Ordinal bytes,
564  char charBuffer[]) {
566  return DSerT::serialize(count, buffer, bytes, charBuffer);
567  else
568  return STI::serialize(count, buffer, bytes, charBuffer);
569  }
570 
573  const char charBuffer[]) {
575  return DSerT::fromIndirectBytesToCount(bytes, charBuffer);
576  else
577  return STI::fromIndirectBytesToCount(bytes, charBuffer);
578  }
579 
581  static void deserialize (const Ordinal bytes,
582  const char charBuffer[],
583  const Ordinal count,
584  VecType buffer[]) {
586  return DSerT::deserialize(bytes, charBuffer, count, buffer);
587  else
588  return STI::deserialize(bytes, charBuffer, count, buffer);
589  }
590 
592 
593  };
594 
596  template <typename Ordinal, typename VecType, typename ValueSerializer>
598 
599  private:
600 
603 
606 
609 
610  public:
611 
613  typedef ValueSerializer value_serializer_type;
614 
616  static const bool supportsDirectSerialization =
618 
621  Ordinal sz_ = 0) :
622  vs(vs_), sz(sz_) {}
623 
625  Ordinal getSerializerSize() const { return sz; }
626 
629  return vs; }
630 
632 
633 
636  const VecType buffer[]) const {
637  return Imp::fromCountToIndirectBytes(*vs, count, buffer, sz);
638  }
639 
641  void serialize (const Ordinal count,
642  const VecType buffer[],
643  const Ordinal bytes,
644  char charBuffer[]) const {
645  Imp::serialize(*vs, count, buffer, bytes, charBuffer, sz);
646  }
647 
650  const char charBuffer[]) const {
651  return Imp::fromIndirectBytesToCount(*vs, bytes, charBuffer, sz);
652  }
653 
655  void deserialize (const Ordinal bytes,
656  const char charBuffer[],
657  const Ordinal count,
658  VecType buffer[]) const {
659  return Imp::deserialize(*vs, bytes, charBuffer, count, buffer, sz);
660  }
661 
663 
664  };
665 
666  }
667 
668 }
669 
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*.
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 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 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)
static Ordinal fromCountToDirectBytes(const Ordinal count)
Return the number of bytes for count objects.
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[])
SerializationImp< Ordinal, VecType, ValueSerializer > Imp
Implementation.
static VecType * convertFromCharPtr(char *ptr)
Convert the pointer type from char*.
Sacado::ValueType< VecType >::type ValueT
Value type.
static ScalarType squareroot(const ScalarType &x)
Teuchos::SerializationTraits< Ordinal, int > iSerT
How to serialize ints.
Ordinal fromCountToIndirectBytes(const Ordinal count, const VecType buffer[]) const
Return the number of bytes for count objects.
DS::DefaultSerializerType ValueSerializer
Default serializer type for values.
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 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.
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 ScalarType conjugate(const ScalarType &x)
Sacado::MP::Vector< storage_double_type > doublePrecision
Serialization implementation for all Vector types.
Ordinal getSerializerSize() const
Return specified serializer size.
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.
Sacado::MP::Vector< storage_half_type > halfPrecision
static Ordinal fromCountToIndirectBytes(const Ordinal count, const VecType buffer[])
Return the number of bytes for count objects.
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.
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.
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
Sacado::mpl::apply< S, ordinal_type, value_double_type >::type storage_double_type
static ScalarType imag(const ScalarType &x)
Teuchos::ScalarTraits< value_type >::coordinateType coordinateType
void deserialize(const Ordinal bytes, const char charBuffer[], const Ordinal count, VecType buffer[]) const
Deserialize from an indirect char[] buffer.
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 const bool supportsDirectSerialization
Whether we support direct serialization.
KOKKOS_INLINE_FUNCTION PCE< Storage > log10(const PCE< Storage > &a)
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)
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.