Eigen  5.0.1-dev+60122df6
 
Loading...
Searching...
No Matches
Serializer.h
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
19namespace Eigen {
20
26template <typename T, typename EnableIf = void>
28
29// Specialization for POD types.
30template <typename T>
31class Serializer<T,
32 typename std::enable_if_t<std::is_trivially_copyable<T>::value && std::is_standard_layout<T>::value>> {
33 public:
40 EIGEN_DEVICE_FUNC size_t size(const T& value) const { return sizeof(value); }
41
49 EIGEN_DEVICE_FUNC uint8_t* serialize(uint8_t* dest, uint8_t* end, const T& value) {
50 if (EIGEN_PREDICT_FALSE(dest == nullptr)) return nullptr;
51 if (EIGEN_PREDICT_FALSE(dest + sizeof(value) > end)) return nullptr;
52 EIGEN_USING_STD(memcpy)
53 memcpy(dest, &value, sizeof(value));
54 return dest + sizeof(value);
55 }
56
64 EIGEN_DEVICE_FUNC const uint8_t* deserialize(const uint8_t* src, const uint8_t* end, T& value) const {
65 if (EIGEN_PREDICT_FALSE(src == nullptr)) return nullptr;
66 if (EIGEN_PREDICT_FALSE(src + sizeof(value) > end)) return nullptr;
67 EIGEN_USING_STD(memcpy)
68 memcpy(&value, src, sizeof(value));
69 return src + sizeof(value);
70 }
71};
72
73// Specialization for DenseBase.
74// Serializes [rows, cols, data...].
75template <typename Derived>
76class Serializer<DenseBase<Derived>, void> {
77 public:
78 typedef typename Derived::Scalar Scalar;
79
80 struct Header {
81 typename Derived::Index rows;
82 typename Derived::Index cols;
83 };
84
85 EIGEN_DEVICE_FUNC size_t size(const Derived& value) const { return sizeof(Header) + sizeof(Scalar) * value.size(); }
86
87 EIGEN_DEVICE_FUNC uint8_t* serialize(uint8_t* dest, uint8_t* end, const Derived& value) {
88 if (EIGEN_PREDICT_FALSE(dest == nullptr)) return nullptr;
89 if (EIGEN_PREDICT_FALSE(dest + size(value) > end)) return nullptr;
90 const size_t header_bytes = sizeof(Header);
91 const size_t data_bytes = sizeof(Scalar) * value.size();
92 Header header = {value.rows(), value.cols()};
93 EIGEN_USING_STD(memcpy)
94 memcpy(dest, &header, header_bytes);
95 dest += header_bytes;
96 memcpy(dest, value.data(), data_bytes);
97 return dest + data_bytes;
98 }
99
100 EIGEN_DEVICE_FUNC const uint8_t* deserialize(const uint8_t* src, const uint8_t* end, Derived& value) const {
101 if (EIGEN_PREDICT_FALSE(src == nullptr)) return nullptr;
102 if (EIGEN_PREDICT_FALSE(src + sizeof(Header) > end)) return nullptr;
103 const size_t header_bytes = sizeof(Header);
104 Header header;
105 EIGEN_USING_STD(memcpy)
106 memcpy(&header, src, header_bytes);
107 src += header_bytes;
108 const size_t data_bytes = sizeof(Scalar) * header.rows * header.cols;
109 if (EIGEN_PREDICT_FALSE(src + data_bytes > end)) return nullptr;
110 value.resize(header.rows, header.cols);
111 memcpy(value.data(), src, data_bytes);
112 return src + data_bytes;
113 }
114};
115
116template <typename Scalar, int Rows, int Cols, int Options, int MaxRows, int MaxCols>
117class Serializer<Matrix<Scalar, Rows, Cols, Options, MaxRows, MaxCols>>
118 : public Serializer<DenseBase<Matrix<Scalar, Rows, Cols, Options, MaxRows, MaxCols>>> {};
119
120template <typename Scalar, int Rows, int Cols, int Options, int MaxRows, int MaxCols>
121class Serializer<Array<Scalar, Rows, Cols, Options, MaxRows, MaxCols>>
122 : public Serializer<DenseBase<Array<Scalar, Rows, Cols, Options, MaxRows, MaxCols>>> {};
123
124namespace internal {
125
126// Recursive serialization implementation helper.
127template <size_t N, typename... Types>
128struct serialize_impl;
129
130template <size_t N, typename T1, typename... Ts>
131struct serialize_impl<N, T1, Ts...> {
132 using Serializer = Eigen::Serializer<typename std::decay<T1>::type>;
133
134 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE size_t serialize_size(const T1& value, const Ts&... args) {
135 Serializer serializer;
136 size_t size = serializer.size(value);
137 return size + serialize_impl<N - 1, Ts...>::serialize_size(args...);
138 }
139
140 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE uint8_t* serialize(uint8_t* dest, uint8_t* end, const T1& value,
141 const Ts&... args) {
142 Serializer serializer;
143 dest = serializer.serialize(dest, end, value);
144 return serialize_impl<N - 1, Ts...>::serialize(dest, end, args...);
145 }
146
147 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const uint8_t* deserialize(const uint8_t* src, const uint8_t* end,
148 T1& value, Ts&... args) {
149 Serializer serializer;
150 src = serializer.deserialize(src, end, value);
151 return serialize_impl<N - 1, Ts...>::deserialize(src, end, args...);
152 }
153};
154
155// Base case.
156template <>
157struct serialize_impl<0> {
158 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE size_t serialize_size() { return 0; }
159
160 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE uint8_t* serialize(uint8_t* dest, uint8_t* /*end*/) { return dest; }
161
162 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const uint8_t* deserialize(const uint8_t* src, const uint8_t* /*end*/) {
163 return src;
164 }
165};
166
167} // namespace internal
168
175template <typename... Args>
176EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE size_t serialize_size(const Args&... args) {
177 return internal::serialize_impl<sizeof...(args), Args...>::serialize_size(args...);
178}
179
188template <typename... Args>
189EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE uint8_t* serialize(uint8_t* dest, uint8_t* end, const Args&... args) {
190 return internal::serialize_impl<sizeof...(args), Args...>::serialize(dest, end, args...);
191}
192
201template <typename... Args>
202EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const uint8_t* deserialize(const uint8_t* src, const uint8_t* end,
203 Args&... args) {
204 return internal::serialize_impl<sizeof...(args), Args...>::deserialize(src, end, args...);
205}
206
207} // namespace Eigen
208
209#endif // EIGEN_SERIALIZER_H
General-purpose arrays with easy API for coefficient-wise operations.
Definition Array.h:48
Base class for all dense matrices, vectors, and arrays.
Definition DenseBase.h:44
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
Definition Serializer.h:27
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
uint8_t * serialize(uint8_t *dest, uint8_t *end, const Args &... args)
Definition Serializer.h:189
size_t serialize_size(const Args &... args)
Definition Serializer.h:176
const uint8_t * deserialize(const uint8_t *src, const uint8_t *end, Args &... args)
Definition Serializer.h:202