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 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef SACADO_MP_SCALAR_TRAITS_IMP_HPP
11 #define SACADO_MP_SCALAR_TRAITS_IMP_HPP
12 
13 #include "Teuchos_ScalarTraits.hpp"
16 #include "Teuchos_Assert.hpp"
18 #include "Sacado_mpl_apply.hpp"
19 
20 namespace Sacado {
21  namespace MP {
22 
23  template <typename S, bool reduct_across_vector>
24  struct ScalarTraitsImp {};
25 
26  // Implementation of Teuchos::ScalarTraits where reductions are taken
27  // across the components of MP::Vector. In this case magnitudeType is
28  // a scalar
29  template <typename S>
30  struct ScalarTraitsImp<S,true> {
32  typedef typename S::value_type value_type;
33  typedef typename S::ordinal_type ordinal_type;
35 
39 
40  typedef typename Sacado::mpl::apply<S,ordinal_type,value_mag_type>::type storage_mag_type;
41  typedef typename Sacado::mpl::apply<S,ordinal_type,value_half_type>::type storage_half_type;
42  typedef typename Sacado::mpl::apply<S,ordinal_type,value_double_type>::type storage_double_type;
43 
48 
49  static const bool isComplex = TVT::isComplex;
50  static const bool isOrdinal = TVT::isOrdinal;
51  static const bool isComparable = TVT::isComparable;
52  static const bool hasMachineParameters = TVT::hasMachineParameters;
53 
54  static value_mag_type eps() { return TVT::eps(); }
55 
56  static value_mag_type sfmin() { return TVT::sfmin(); }
57 
58  static value_mag_type base() { return TVT::base(); }
59 
60  static value_mag_type prec() { return TVT::prec(); }
61 
62  static value_mag_type t() { return TVT::t(); }
63 
64  static value_mag_type rnd() { return TVT::rnd(); }
65 
66  static value_mag_type emin() { return TVT::emin(); }
67 
68  static value_mag_type rmin() { return TVT::rmin(); }
69 
70  static value_mag_type emax() { return TVT::emax(); }
71 
72  static value_mag_type rmax() { return TVT::rmax(); }
73 
74  static magnitudeType magnitude(const ScalarType& a) {
76  const ordinal_type sz = a.size();
77  for (ordinal_type i=0; i<sz; ++i) {
78  value_mag_type t = TVT::magnitude(a.fastAccessCoeff(i));
79  m +=t*t;
80  }
81  return std::sqrt(m);
82  }
83 
84  static ScalarType zero() { return ScalarType(0.0); }
85 
86  static ScalarType one() { return ScalarType(1.0); }
87 
88 
89  static ScalarType conjugate(const ScalarType& x) {
90  int sz = x.size();
91  ScalarType y(sz, value_type(0.0));
92  for (int i=0; i<sz; i++)
93  y.fastAccessCoeff(i) = TVT::conjugate(x.fastAccessCoeff(i));
94  return y;
95  }
96 
97 
98  static magnitudeType real(const ScalarType& x) {
100  const ordinal_type sz = x.size();
101  for (ordinal_type i=0; i<sz; ++i) {
102  value_mag_type t = TVT::real(x.fastAccessCoeff(i));
103  m +=t*t;
104  }
105  return std::sqrt(m);
106  }
107 
108 
109  static magnitudeType imag(const ScalarType& x) {
110  magnitudeType m = magnitudeType(0.0);
111  const ordinal_type sz = x.size();
112  for (ordinal_type i=0; i<sz; ++i) {
113  value_mag_type t = TVT::imag(x.fastAccessCoeff(i));
114  m +=t*t;
115  }
116  return std::sqrt(m);
117  }
118 
119  static value_type nan() { return TVT::nan(); }
120 
121  static bool isnaninf(const ScalarType& x) {
122  for (int i=0; i<x.size(); i++)
123  if (TVT::isnaninf(x.fastAccessCoeff(i)))
124  return true;
125  return false;
126  }
127 
128  static void seedrandom(unsigned int s) { TVT::seedrandom(s); }
129 
130  static ScalarType random() { return ScalarType(TVT::random()); }
131 
132  static const char * name() { return "Sacado::MP::Vector<>"; }
133 
134  static ScalarType squareroot(const ScalarType& x) { return std::sqrt(x); }
135 
136  static ScalarType pow(const ScalarType& x, const ScalarType& y) {
137  return std::pow(x,y);
138  }
139 
140  static ScalarType log(const ScalarType& x) { return std::log(x); }
141 
142  static ScalarType log10(const ScalarType& x) { return std::log10(x); }
143 
144  }; // class ScalarTraitsImp<S,true>
145 
146  // Implementation of Teuchos::ScalarTraits where reductions are not taken
147  // across the components of MP::Vector. In this case magnitudeType is
148  // an MP::Vector
149  template <typename S>
150  struct ScalarTraitsImp<S,false> {
152  typedef typename S::value_type value_type;
153  typedef typename S::ordinal_type ordinal_type;
155 
159 
160  typedef typename Sacado::mpl::apply<S,ordinal_type,value_mag_type>::type storage_mag_type;
161  typedef typename Sacado::mpl::apply<S,ordinal_type,value_half_type>::type storage_half_type;
162  typedef typename Sacado::mpl::apply<S,ordinal_type,value_double_type>::type storage_double_type;
163 
168 
169  static const bool isComplex = TVT::isComplex;
170  static const bool isOrdinal = TVT::isOrdinal;
171  static const bool isComparable = TVT::isComparable;
172  static const bool hasMachineParameters = TVT::hasMachineParameters;
173 
174  static value_mag_type eps() { return TVT::eps(); }
175 
176  static value_mag_type sfmin() { return TVT::sfmin(); }
177 
178  static value_mag_type base() { return TVT::base(); }
179 
180  static value_mag_type prec() { return TVT::prec(); }
181 
182  static value_mag_type t() { return TVT::t(); }
183 
184  static value_mag_type rnd() { return TVT::rnd(); }
185 
186  static value_mag_type emin() { return TVT::emin(); }
187 
188  static value_mag_type rmin() { return TVT::rmin(); }
189 
190  static value_mag_type emax() { return TVT::emax(); }
191 
192  static value_mag_type rmax() { return TVT::rmax(); }
193 
194  static magnitudeType magnitude(const ScalarType& a) {
195  return std::fabs(a);
196  }
197 
198  static ScalarType zero() { return ScalarType(0.0); }
199 
200  static ScalarType one() { return ScalarType(1.0); }
201 
202 
203  static ScalarType conjugate(const ScalarType& x) {
204  int sz = x.size();
205  ScalarType y(sz, value_type(0.0));
206  for (int i=0; i<sz; i++)
207  y.fastAccessCoeff(i) = TVT::conjugate(x.fastAccessCoeff(i));
208  return y;
209  }
210 
211  static ScalarType real(const ScalarType& x) {
212  int sz = x.size();
213  ScalarType y(sz, value_type(0.0));
214  for (int i=0; i<sz; i++)
215  y.fastAccessCoeff(i) = TVT::real(x.fastAccessCoeff(i));
216  return y;
217  }
218 
219  static ScalarType imag(const ScalarType& x) {
220  int sz = x.size();
221  ScalarType y(sz, value_type(0.0));
222  for (int i=0; i<sz; i++)
223  y.fastAccessCoeff(i) = TVT::imag(x.fastAccessCoeff(i));
224  return y;
225  }
226 
227  static value_type nan() { return TVT::nan(); }
228 
229  static bool isnaninf(const ScalarType& x) {
230  for (int i=0; i<x.size(); i++)
231  if (TVT::isnaninf(x.fastAccessCoeff(i)))
232  return true;
233  return false;
234  }
235 
236  static void seedrandom(unsigned int s) { TVT::seedrandom(s); }
237 
238  static ScalarType random() { return ScalarType(TVT::random()); }
239 
240  static const char * name() { return "Sacado::MP::Vector<>"; }
241 
242  static ScalarType squareroot(const ScalarType& x) { return std::sqrt(x); }
243 
244  static ScalarType pow(const ScalarType& x, const ScalarType& y) {
245  return std::pow(x,y);
246  }
247 
248  static ScalarType log(const ScalarType& x) { return std::log(x); }
249 
250  static ScalarType log10(const ScalarType& x) { return std::log10(x); }
251 
252  }; // class ScalarTraitsImp<S,false>
253 
255  template <typename Ordinal, typename VecType, typename Serializer>
257 
258  private:
259 
261  typedef typename Sacado::ValueType<VecType>::type ValueT;
262 
265 
268 
269  public:
270 
272  static const bool supportsDirectSerialization = false;
273 
275 
276 
278  static Ordinal fromCountToIndirectBytes(const Serializer& vs,
279  const Ordinal count,
280  const VecType buffer[],
281  const Ordinal sz = 0) {
282  Ordinal bytes = 0;
283  VecType *x = NULL;
284  const VecType *cx;
285  for (Ordinal i=0; i<count; i++) {
286  int my_sz = buffer[i].size();
287  int tot_sz = sz;
288  if (sz == 0) tot_sz = my_sz;
289  Ordinal b1 = iSerT::fromCountToIndirectBytes(1, &tot_sz);
290  if (tot_sz != my_sz) {
291  if (x == NULL)
292  x = new VecType;
293  *x = buffer[i];
294  x->reset(tot_sz);
295  cx = x;
296  }
297  else
298  cx = &(buffer[i]);
299  Ordinal b2 = vs.fromCountToIndirectBytes(tot_sz, cx->coeff());
301  bytes += b1+b2+b3;
302  }
303  if (x != NULL)
304  delete x;
305  return bytes;
306  }
307 
309  static void serialize (const Serializer& vs,
310  const Ordinal count,
311  const VecType buffer[],
312  const Ordinal bytes,
313  char charBuffer[],
314  const Ordinal sz = 0) {
315  VecType *x = NULL;
316  const VecType *cx;
317  for (Ordinal i=0; i<count; i++) {
318  // First serialize size
319  int my_sz = buffer[i].size();
320  int tot_sz = sz;
321  if (sz == 0) tot_sz = my_sz;
322  Ordinal b1 = iSerT::fromCountToIndirectBytes(1, &tot_sz);
323  iSerT::serialize(1, &tot_sz, b1, charBuffer);
324  charBuffer += b1;
325 
326  // Next serialize vector coefficients
327  if (tot_sz != my_sz) {
328  if (x == NULL)
329  x = new VecType;
330  *x = buffer[i];
331  x->reset(tot_sz);
332  cx = x;
333  }
334  else
335  cx = &(buffer[i]);
336  Ordinal b2 = vs.fromCountToIndirectBytes(tot_sz, cx->coeff());
338  oSerT::serialize(1, &b2, b3, charBuffer);
339  charBuffer += b3;
340  vs.serialize(tot_sz, cx->coeff(), b2, charBuffer);
341  charBuffer += b2;
342  }
343  if (x != NULL)
344  delete x;
345  }
346 
348  static Ordinal fromIndirectBytesToCount(const Serializer& vs,
349  const Ordinal bytes,
350  const char charBuffer[],
351  const Ordinal sz = 0) {
352  Ordinal count = 0;
353  Ordinal bytes_used = 0;
354  while (bytes_used < bytes) {
355 
356  // Bytes for size
358  bytes_used += b1;
359  charBuffer += b1;
360 
361  // Bytes for vector coefficients
363  const Ordinal *b2 = oSerT::convertFromCharPtr(charBuffer);
364  bytes_used += b3;
365  charBuffer += b3;
366  bytes_used += *b2;
367  charBuffer += *b2;
368 
369  ++count;
370  }
371  return count;
372  }
373 
375  static void deserialize (const Serializer& vs,
376  const Ordinal bytes,
377  const char charBuffer[],
378  const Ordinal count,
379  VecType buffer[],
380  const Ordinal sz = 0) {
381  for (Ordinal i=0; i<count; i++) {
382 
383  // Deserialize size
385  const int *my_sz = iSerT::convertFromCharPtr(charBuffer);
386  charBuffer += b1;
387 
388  // Create empty Vector object of given size
389  int tot_sz = sz;
390  if (sz == 0) tot_sz = *my_sz;
391  buffer[i] = VecType(tot_sz, ValueT(0.0));
392 
393  // Deserialize vector coefficients
395  const Ordinal *b2 = oSerT::convertFromCharPtr(charBuffer);
396  charBuffer += b3;
397  vs.deserialize(*b2, charBuffer, *my_sz, buffer[i].coeff());
398  charBuffer += *b2;
399  }
400 
401  }
403 
404  };
405 
407  template <typename Ordinal, typename VecType, bool is_static = false>
409 
410  private:
411 
413  typedef typename Sacado::ValueType<VecType>::type ValueT;
414 
417 
420 
423 
424  public:
425 
427  static const bool supportsDirectSerialization =
429 
431 
432 
435  const VecType buffer[]) {
437  DS::getDefaultSerializer(), count, buffer);
438  }
439 
441  static void serialize (const Ordinal count,
442  const VecType buffer[],
443  const Ordinal bytes,
444  char charBuffer[]) {
446  DS::getDefaultSerializer(), count, buffer, bytes, charBuffer);
447  }
448 
451  const char charBuffer[]) {
453  DS::getDefaultSerializer(), bytes, charBuffer);
454  }
455 
457  static void deserialize (const Ordinal bytes,
458  const char charBuffer[],
459  const Ordinal count,
460  VecType buffer[]) {
462  DS::getDefaultSerializer(), bytes, charBuffer, count, buffer);
463  }
464 
466 
467  };
468 
470  template <typename Ordinal, typename VecType>
471  struct SerializationTraitsImp<Ordinal, VecType, true> {
472  typedef typename Sacado::ValueType<VecType>::type ValueT;
476 
478  static const bool supportsDirectSerialization =
479  vSerT::supportsDirectSerialization;
480 
482 
483 
485  static Ordinal fromCountToDirectBytes(const Ordinal count) {
486  return DSerT::fromCountToDirectBytes(count);
487  }
488 
490  static char* convertToCharPtr( VecType* ptr ) {
491  return DSerT::convertToCharPtr(ptr);
492  }
493 
495  static const char* convertToCharPtr( const VecType* ptr ) {
496  return DSerT::convertToCharPtr(ptr);
497  }
498 
500  static Ordinal fromDirectBytesToCount(const Ordinal bytes) {
501  return DSerT::fromDirectBytesToCount(bytes);
502  }
503 
505  static VecType* convertFromCharPtr( char* ptr ) {
506  return DSerT::convertFromCharPtr(ptr);
507  }
508 
510  static const VecType* convertFromCharPtr( const char* ptr ) {
511  return DSerT::convertFromCharPtr(ptr);
512  }
513 
515 
517 
518 
521  const VecType buffer[]) {
523  return DSerT::fromCountToIndirectBytes(count, buffer);
524  else
525  return STI::fromCountToIndirectBytes(count, buffer);
526  }
527 
529  static void serialize (const Ordinal count,
530  const VecType buffer[],
531  const Ordinal bytes,
532  char charBuffer[]) {
534  return DSerT::serialize(count, buffer, bytes, charBuffer);
535  else
536  return STI::serialize(count, buffer, bytes, charBuffer);
537  }
538 
541  const char charBuffer[]) {
543  return DSerT::fromIndirectBytesToCount(bytes, charBuffer);
544  else
545  return STI::fromIndirectBytesToCount(bytes, charBuffer);
546  }
547 
549  static void deserialize (const Ordinal bytes,
550  const char charBuffer[],
551  const Ordinal count,
552  VecType buffer[]) {
554  return DSerT::deserialize(bytes, charBuffer, count, buffer);
555  else
556  return STI::deserialize(bytes, charBuffer, count, buffer);
557  }
558 
560 
561  };
562 
564  template <typename Ordinal, typename VecType, typename ValueSerializer>
566 
567  private:
568 
571 
574 
577 
578  public:
579 
581  typedef ValueSerializer value_serializer_type;
582 
584  static const bool supportsDirectSerialization =
586 
589  Ordinal sz_ = 0) :
590  vs(vs_), sz(sz_) {}
591 
593  Ordinal getSerializerSize() const { return sz; }
594 
597  return vs; }
598 
600 
601 
604  const VecType buffer[]) const {
605  return Imp::fromCountToIndirectBytes(*vs, count, buffer, sz);
606  }
607 
609  void serialize (const Ordinal count,
610  const VecType buffer[],
611  const Ordinal bytes,
612  char charBuffer[]) const {
613  Imp::serialize(*vs, count, buffer, bytes, charBuffer, sz);
614  }
615 
618  const char charBuffer[]) const {
619  return Imp::fromIndirectBytesToCount(*vs, bytes, charBuffer, sz);
620  }
621 
623  void deserialize (const Ordinal bytes,
624  const char charBuffer[],
625  const Ordinal count,
626  VecType buffer[]) const {
627  return Imp::deserialize(*vs, bytes, charBuffer, count, buffer, sz);
628  }
629 
631 
632  };
633 
634  }
635 
636 }
637 
638 #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.