Sacado Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Sacado_Fad_SerializationTraitsImp.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Sacado Package
4 //
5 // Copyright 2006 NTESS and the Sacado contributors.
6 // SPDX-License-Identifier: LGPL-2.1-or-later
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef SACADO_FAD_SERIALIZATIONTRAITSIMP_HPP
11 #define SACADO_FAD_SERIALIZATIONTRAITSIMP_HPP
12 
13 #include "Sacado_ConfigDefs.h"
14 
15 #ifdef HAVE_SACADO_TEUCHOSCOMM
16 
19 #include "Teuchos_RCP.hpp"
20 
21 namespace Sacado {
22 
23  namespace Fad {
24 
26  template <typename Ordinal, typename FadType, typename Serializer>
27  struct SerializationImp {
28 
29  private:
30 
33 
36 
38  typedef typename Sacado::ValueType<FadType>::type value_type;
39 
40  public:
41 
43  static const bool supportsDirectSerialization = false;
44 
46 
47 
49  static Ordinal fromCountToIndirectBytes(const Serializer& vs,
50  const Ordinal count,
51  const FadType buffer[],
52  const Ordinal sz = 0) {
53  Ordinal bytes = 0;
54  FadType *x = NULL;
55  const FadType *cx;
56  for (Ordinal i=0; i<count; i++) {
57  int my_sz = buffer[i].size();
58  int tot_sz = sz;
59  if (sz == 0) tot_sz = my_sz;
60  Ordinal b1 = iSerT::fromCountToIndirectBytes(1, &tot_sz);
61  Ordinal b2 = vs.fromCountToIndirectBytes(1, &(buffer[i].val()));
62  Ordinal b3 = oSerT::fromCountToIndirectBytes(1, &b2);
63  Ordinal b4;
64  if (tot_sz != my_sz) {
65  if (x == NULL)
66  x = new FadType(tot_sz, 0.0);
67  *x = buffer[i];
68  x->expand(tot_sz);
69  cx = x;
70  }
71  else
72  cx = &(buffer[i]);
73  b4 = vs.fromCountToIndirectBytes(tot_sz, cx->dx());
74  Ordinal b5 = oSerT::fromCountToIndirectBytes(1, &b4);
75  bytes += b1+b2+b3+b4+b5;
76  }
77  if (x != NULL)
78  delete x;
79  return bytes;
80  }
81 
83  static void serialize (const Serializer& vs,
84  const Ordinal count,
85  const FadType buffer[],
86  const Ordinal bytes,
87  char charBuffer[],
88  const Ordinal sz = 0) {
89  FadType *x = NULL;
90  const FadType *cx;
91  for (Ordinal i=0; i<count; i++) {
92  // First serialize size
93  int my_sz = buffer[i].size();
94  int tot_sz = sz;
95  if (sz == 0) tot_sz = my_sz;
96  Ordinal b1 = iSerT::fromCountToIndirectBytes(1, &tot_sz);
97  iSerT::serialize(1, &tot_sz, b1, charBuffer);
98  charBuffer += b1;
99 
100  // Next serialize value
101  Ordinal b2 = vs.fromCountToIndirectBytes(1, &(buffer[i].val()));
102  Ordinal b3 = oSerT::fromCountToIndirectBytes(1, &b2);
103  oSerT::serialize(1, &b2, b3, charBuffer);
104  charBuffer += b3;
105  vs.serialize(1, &(buffer[i].val()), b2, charBuffer);
106  charBuffer += b2;
107 
108  // Next serialize derivative array
109  Ordinal b4;
110  if (tot_sz != my_sz) {
111  if (x == NULL)
112  x = new FadType(tot_sz, 0.0);
113  *x = buffer[i];
114  x->expand(tot_sz);
115  cx = x;
116  }
117  else
118  cx = &(buffer[i]);
119  b4 = vs.fromCountToIndirectBytes(tot_sz, cx->dx());
120  Ordinal b5 = oSerT::fromCountToIndirectBytes(1, &b4);
121  oSerT::serialize(1, &b4, b5, charBuffer);
122  charBuffer += b5;
123  vs.serialize(tot_sz, cx->dx(), b4, charBuffer);
124  charBuffer += b4;
125  }
126  if (x != NULL)
127  delete x;
128  }
129 
131  static Ordinal fromIndirectBytesToCount(const Serializer& vs,
132  const Ordinal bytes,
133  const char charBuffer[],
134  const Ordinal sz = 0) {
135  Ordinal count = 0;
136  Ordinal bytes_used = 0;
137  while (bytes_used < bytes) {
138 
139  // Bytes for size
140  Ordinal b1 = iSerT::fromCountToDirectBytes(1);
141  bytes_used += b1;
142  charBuffer += b1;
143 
144  // Bytes for value
145  Ordinal b3 = oSerT::fromCountToDirectBytes(1);
146  const Ordinal *b2 = oSerT::convertFromCharPtr(charBuffer);
147  bytes_used += b3;
148  charBuffer += b3;
149  bytes_used += *b2;
150  charBuffer += *b2;
151 
152  // Bytes for derivative array
153  Ordinal b5 = oSerT::fromCountToDirectBytes(1);
154  const Ordinal *b4 = oSerT::convertFromCharPtr(charBuffer);
155  bytes_used += b5;
156  charBuffer += b5;
157  bytes_used += *b4;
158  charBuffer += *b4;
159 
160  ++count;
161  }
162  return count;
163  }
164 
166  static void deserialize (const Serializer& vs,
167  const Ordinal bytes,
168  const char charBuffer[],
169  const Ordinal count,
170  FadType buffer[],
171  const Ordinal sz = 0) {
172  for (Ordinal i=0; i<count; i++) {
173 
174  // Deserialize size
175  Ordinal b1 = iSerT::fromCountToDirectBytes(1);
176  const int *my_sz = iSerT::convertFromCharPtr(charBuffer);
177  charBuffer += b1;
178 
179  // Create empty Fad object of given size
180  int tot_sz = sz;
181  if (sz == 0) tot_sz = *my_sz;
182  buffer[i] = FadType(tot_sz, 0.0);
183 
184  // Deserialize value
185  Ordinal b3 = oSerT::fromCountToDirectBytes(1);
186  const Ordinal *b2 = oSerT::convertFromCharPtr(charBuffer);
187  charBuffer += b3;
188  vs.deserialize(*b2, charBuffer, 1, &(buffer[i].val()));
189  charBuffer += *b2;
190 
191  // Deserialize derivative array
192  Ordinal b5 = oSerT::fromCountToDirectBytes(1);
193  const Ordinal *b4 = oSerT::convertFromCharPtr(charBuffer);
194  charBuffer += b5;
195  vs.deserialize(*b4, charBuffer, *my_sz,
196  &(buffer[i].fastAccessDx(0)));
197  charBuffer += *b4;
198  }
199 
200  }
201 
203 
204  };
205 
207  template <typename Ordinal, typename FadType>
208  struct SerializationTraitsImp {
209 
210  private:
211 
213  typedef typename Sacado::ValueType<FadType>::type ValueT;
214 
217 
219  typedef typename DS::DefaultSerializerType ValueSerializer;
220 
222  typedef SerializationImp<Ordinal,FadType,ValueSerializer> Imp;
223 
224  public:
225 
227  static const bool supportsDirectSerialization =
228  Imp::supportsDirectSerialization;
229 
231 
232 
234  static Ordinal fromCountToIndirectBytes(const Ordinal count,
235  const FadType buffer[]) {
236  return Imp::fromCountToIndirectBytes(
237  DS::getDefaultSerializer(), count, buffer);
238  }
239 
241  static void serialize (const Ordinal count,
242  const FadType buffer[],
243  const Ordinal bytes,
244  char charBuffer[]) {
245  Imp::serialize(
246  DS::getDefaultSerializer(), count, buffer, bytes, charBuffer);
247  }
248 
250  static Ordinal fromIndirectBytesToCount(const Ordinal bytes,
251  const char charBuffer[]) {
252  return Imp::fromIndirectBytesToCount(
253  DS::getDefaultSerializer(), bytes, charBuffer);
254  }
255 
257  static void deserialize (const Ordinal bytes,
258  const char charBuffer[],
259  const Ordinal count,
260  FadType buffer[]) {
261  Imp::deserialize(
262  DS::getDefaultSerializer(), bytes, charBuffer, count, buffer);
263  }
264 
266 
267  };
268 
270  template <typename Ordinal, typename FadType>
271  struct StaticSerializationTraitsImp {
272  typedef typename Sacado::ValueType<FadType>::type ValueT;
275  typedef Sacado::Fad::SerializationTraitsImp<Ordinal,FadType> STI;
276 
278  static const bool supportsDirectSerialization =
279  vSerT::supportsDirectSerialization;
280 
282 
283 
285  static Ordinal fromCountToDirectBytes(const Ordinal count) {
286  return DSerT::fromCountToDirectBytes(count);
287  }
288 
290  static char* convertToCharPtr( FadType* ptr ) {
291  return DSerT::convertToCharPtr(ptr);
292  }
293 
295  static const char* convertToCharPtr( const FadType* ptr ) {
296  return DSerT::convertToCharPtr(ptr);
297  }
298 
300  static Ordinal fromDirectBytesToCount(const Ordinal bytes) {
301  return DSerT::fromDirectBytesToCount(bytes);
302  }
303 
305  static FadType* convertFromCharPtr( char* ptr ) {
306  return DSerT::convertFromCharPtr(ptr);
307  }
308 
310  static const FadType* convertFromCharPtr( const char* ptr ) {
311  return DSerT::convertFromCharPtr(ptr);
312  }
313 
315 
317 
318 
320  static Ordinal fromCountToIndirectBytes(const Ordinal count,
321  const FadType buffer[]) {
322  if (supportsDirectSerialization)
323  return DSerT::fromCountToIndirectBytes(count, buffer);
324  else
325  return STI::fromCountToIndirectBytes(count, buffer);
326  }
327 
329  static void serialize (const Ordinal count,
330  const FadType buffer[],
331  const Ordinal bytes,
332  char charBuffer[]) {
333  if (supportsDirectSerialization)
334  return DSerT::serialize(count, buffer, bytes, charBuffer);
335  else
336  return STI::serialize(count, buffer, bytes, charBuffer);
337  }
338 
340  static Ordinal fromIndirectBytesToCount(const Ordinal bytes,
341  const char charBuffer[]) {
342  if (supportsDirectSerialization)
343  return DSerT::fromIndirectBytesToCount(bytes, charBuffer);
344  else
345  return STI::fromIndirectBytesToCount(bytes, charBuffer);
346  }
347 
349  static void deserialize (const Ordinal bytes,
350  const char charBuffer[],
351  const Ordinal count,
352  FadType buffer[]) {
353  if (supportsDirectSerialization)
354  return DSerT::deserialize(bytes, charBuffer, count, buffer);
355  else
356  return STI::deserialize(bytes, charBuffer, count, buffer);
357  }
358 
360 
361  };
362 
364  template <typename Ordinal, typename FadType, typename ValueSerializer>
365  class SerializerImp {
366 
367  private:
368 
370  typedef SerializationImp<Ordinal,FadType,ValueSerializer> Imp;
371 
374 
376  Ordinal sz;
377 
378  public:
379 
381  typedef ValueSerializer value_serializer_type;
382 
384  static const bool supportsDirectSerialization =
385  Imp::supportsDirectSerialization;
386 
388  SerializerImp(const Teuchos::RCP<const ValueSerializer>& vs_,
389  Ordinal sz_ = 0) :
390  vs(vs_), sz(sz_) {}
391 
393  Ordinal getSerializerSize() const { return sz; }
394 
396  Teuchos::RCP<const value_serializer_type> getValueSerializer() const {
397  return vs; }
398 
400 
401 
403  Ordinal fromCountToIndirectBytes(const Ordinal count,
404  const FadType buffer[]) const {
405  return Imp::fromCountToIndirectBytes(*vs, count, buffer, sz);
406  }
407 
409  void serialize (const Ordinal count,
410  const FadType buffer[],
411  const Ordinal bytes,
412  char charBuffer[]) const {
413  Imp::serialize(*vs, count, buffer, bytes, charBuffer, sz);
414  }
415 
417  Ordinal fromIndirectBytesToCount(const Ordinal bytes,
418  const char charBuffer[]) const {
419  return Imp::fromIndirectBytesToCount(*vs, bytes, charBuffer, sz);
420  }
421 
423  void deserialize (const Ordinal bytes,
424  const char charBuffer[],
425  const Ordinal count,
426  FadType buffer[]) const {
427  return Imp::deserialize(*vs, bytes, charBuffer, count, buffer, sz);
428  }
429 
431 
432  };
433 
434  } // namespace Fad
435 
436 } // namespace Sacado
437 
438 #endif // HAVE_SACADO_TEUCHOSCOMM
439 
440 #endif // SACADO_FAD_SERIALIZATIONTRAITSIMP_HPP
int * count
Sacado::Fad::DFad< double > FadType
expr val()
int Ordinal
expr expr expr fastAccessDx(i)) FAD_UNARYOP_MACRO(exp