Serializer.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2021 The Eigen Team
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_SERIALIZER_H
11 #define EIGEN_SERIALIZER_H
12 
13 #include <type_traits>
14 
15 // The Serializer class encodes data into a memory buffer so it can be later
16 // reconstructed. This is mainly used to send objects back-and-forth between
17 // the CPU and GPU.
18 
19 namespace Eigen {
20 
26 template<typename T, typename EnableIf = void>
27 class Serializer;
28 
29 // Specialization for POD types.
30 template<typename T>
31 class Serializer<T, typename std::enable_if_t<
32  std::is_trivial<T>::value
33  && std::is_standard_layout<T>::value>> {
34  public:
35 
42  EIGEN_DEVICE_FUNC size_t size(const T& value) const {
43  return sizeof(value);
44  }
45 
54  if (EIGEN_PREDICT_FALSE(dest == nullptr)) return nullptr;
55  if (EIGEN_PREDICT_FALSE(dest + sizeof(value) > end)) return nullptr;
56  EIGEN_USING_STD(memcpy)
57  memcpy(dest, &value, sizeof(value));
58  return dest + sizeof(value);
59  }
60 
68  EIGEN_DEVICE_FUNC const uint8_t* deserialize(const uint8_t* src, const uint8_t* end, T& value) const {
69  if (EIGEN_PREDICT_FALSE(src == nullptr)) return nullptr;
70  if (EIGEN_PREDICT_FALSE(src + sizeof(value) > end)) return nullptr;
71  EIGEN_USING_STD(memcpy)
72  memcpy(&value, src, sizeof(value));
73  return src + sizeof(value);
74  }
75 };
76 
77 // Specialization for DenseBase.
78 // Serializes [rows, cols, data...].
79 template<typename Derived>
80 class Serializer<DenseBase<Derived>, void> {
81  public:
82  typedef typename Derived::Scalar Scalar;
83 
84  struct Header {
85  typename Derived::Index rows;
86  typename Derived::Index cols;
87  };
88 
89  EIGEN_DEVICE_FUNC size_t size(const Derived& value) const {
90  return sizeof(Header) + sizeof(Scalar) * value.size();
91  }
92 
93  EIGEN_DEVICE_FUNC uint8_t* serialize(uint8_t* dest, uint8_t* end, const Derived& value) {
94  if (EIGEN_PREDICT_FALSE(dest == nullptr)) return nullptr;
95  if (EIGEN_PREDICT_FALSE(dest + size(value) > end)) return nullptr;
96  const size_t header_bytes = sizeof(Header);
97  const size_t data_bytes = sizeof(Scalar) * value.size();
98  Header header = {value.rows(), value.cols()};
99  EIGEN_USING_STD(memcpy)
100  memcpy(dest, &header, header_bytes);
101  dest += header_bytes;
102  memcpy(dest, value.data(), data_bytes);
103  return dest + data_bytes;
104  }
105 
106  EIGEN_DEVICE_FUNC const uint8_t* deserialize(const uint8_t* src, const uint8_t* end, Derived& value) const {
107  if (EIGEN_PREDICT_FALSE(src == nullptr)) return nullptr;
108  if (EIGEN_PREDICT_FALSE(src + sizeof(Header) > end)) return nullptr;
109  const size_t header_bytes = sizeof(Header);
110  Header header;
111  EIGEN_USING_STD(memcpy)
112  memcpy(&header, src, header_bytes);
113  src += header_bytes;
114  const size_t data_bytes = sizeof(Scalar) * header.rows * header.cols;
115  if (EIGEN_PREDICT_FALSE(src + data_bytes > end)) return nullptr;
116  value.resize(header.rows, header.cols);
117  memcpy(value.data(), src, data_bytes);
118  return src + data_bytes;
119  }
120 };
121 
122 template<typename Scalar, int Rows, int Cols, int Options, int MaxRows, int MaxCols>
123 class Serializer<Matrix<Scalar, Rows, Cols, Options, MaxRows, MaxCols> > : public
124  Serializer<DenseBase<Matrix<Scalar, Rows, Cols, Options, MaxRows, MaxCols> > > {};
125 
126 template<typename Scalar, int Rows, int Cols, int Options, int MaxRows, int MaxCols>
127 class Serializer<Array<Scalar, Rows, Cols, Options, MaxRows, MaxCols> > : public
128  Serializer<DenseBase<Array<Scalar, Rows, Cols, Options, MaxRows, MaxCols> > > {};
129 
130 namespace internal {
131 
132 // Recursive serialization implementation helper.
133 template<size_t N, typename... Types>
134 struct serialize_impl;
135 
136 template<size_t N, typename T1, typename... Ts>
137 struct serialize_impl<N, T1, Ts...> {
139 
140  static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
141  size_t serialize_size(const T1& value, const Ts&... args) {
142  Serializer serializer;
143  size_t size = serializer.size(value);
145  }
146 
147  static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
148  uint8_t* serialize(uint8_t* dest, uint8_t* end, const T1& value, const Ts&... args) {
149  Serializer serializer;
150  dest = serializer.serialize(dest, end, value);
151  return serialize_impl<N-1, Ts...>::serialize(dest, end, args...);
152  }
153 
154  static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
155  const uint8_t* deserialize(const uint8_t* src, const uint8_t* end, T1& value, Ts&... args) {
156  Serializer serializer;
157  src = serializer.deserialize(src, end, value);
158  return serialize_impl<N-1, Ts...>::deserialize(src, end, args...);
159  }
160 };
161 
162 // Base case.
163 template<>
164 struct serialize_impl<0> {
165  static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
166  size_t serialize_size() { return 0; }
167 
168  static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
169  uint8_t* serialize(uint8_t* dest, uint8_t* /*end*/) { return dest; }
170 
171  static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
172  const uint8_t* deserialize(const uint8_t* src, const uint8_t* /*end*/) { return src; }
173 };
174 
175 } // namespace internal
176 
177 
184 template<typename... Args>
185 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
186 size_t serialize_size(const Args&... args) {
187  return internal::serialize_impl<sizeof...(args), Args...>::serialize_size(args...);
188 }
189 
198 template<typename... Args>
199 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
200 uint8_t* serialize(uint8_t* dest, uint8_t* end, const Args&... args) {
201  return internal::serialize_impl<sizeof...(args), Args...>::serialize(dest, end, args...);
202 }
203 
212 template<typename... Args>
213 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
214 const uint8_t* deserialize(const uint8_t* src, const uint8_t* end, Args&... args) {
215  return internal::serialize_impl<sizeof...(args), Args...>::deserialize(src, end, args...);
216 }
217 
218 } // namespace Eigen
219 
220 #endif // EIGEN_SERIALIZER_H
if((m *x).isApprox(y))
#define EIGEN_USING_STD(FUNC)
Definition: Macros.h:1080
#define EIGEN_PREDICT_FALSE(x)
Definition: Macros.h:1177
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:883
General-purpose arrays with easy API for coefficient-wise operations.
Definition: Array.h:49
Base class for all dense matrices, vectors, and arrays.
Definition: DenseBase.h:42
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182
const uint8_t * deserialize(const uint8_t *src, const uint8_t *end, Derived &value) const
Definition: Serializer.h:106
uint8_t * serialize(uint8_t *dest, uint8_t *end, const Derived &value)
Definition: Serializer.h:93
size_t size(const Derived &value) const
Definition: Serializer.h:89
static const lastp1_t end
std::uint8_t uint8_t
Definition: Meta.h:35
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
uint8_t * serialize(uint8_t *dest, uint8_t *end, const Args &... args)
Definition: Serializer.h:200
size_t serialize_size(const Args &... args)
Definition: Serializer.h:186
const uint8_t * deserialize(const uint8_t *src, const uint8_t *end, Args &... args)
Definition: Serializer.h:214
Definition: BFloat16.h:222