mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-07-18 10:54:26 +08:00
Add a simple serialization mechanism.
The `Serializer<T>` class implements a binary serialization that can write to (`serialize`) and read from (`deserialize`) a byte buffer. Also added convenience routines for serializing a list of arguments. This will mainly be for testing, specifically to transfer data to and from the GPU.
This commit is contained in:
parent
558b3d4fb8
commit
fcd73b4884
@ -166,6 +166,7 @@ using std::ptrdiff_t;
|
|||||||
#include "src/Core/util/XprHelper.h"
|
#include "src/Core/util/XprHelper.h"
|
||||||
#include "src/Core/util/Memory.h"
|
#include "src/Core/util/Memory.h"
|
||||||
#include "src/Core/util/IntegralConstant.h"
|
#include "src/Core/util/IntegralConstant.h"
|
||||||
|
#include "src/Core/util/Serializer.h"
|
||||||
#include "src/Core/util/SymbolicIndex.h"
|
#include "src/Core/util/SymbolicIndex.h"
|
||||||
|
|
||||||
#include "src/Core/NumTraits.h"
|
#include "src/Core/NumTraits.h"
|
||||||
|
207
Eigen/src/Core/util/Serializer.h
Normal file
207
Eigen/src/Core/util/Serializer.h
Normal file
@ -0,0 +1,207 @@
|
|||||||
|
// This file is part of Eigen, a lightweight C++ template library
|
||||||
|
// for linear algebra.
|
||||||
|
//
|
||||||
|
// Copyright (C) 2021 The Eigen Team
|
||||||
|
//
|
||||||
|
// This Source Code Form is subject to the terms of the Mozilla
|
||||||
|
// Public License v. 2.0. If a copy of the MPL was not distributed
|
||||||
|
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
||||||
|
|
||||||
|
#ifndef EIGEN_SERIALIZER_H
|
||||||
|
#define EIGEN_SERIALIZER_H
|
||||||
|
|
||||||
|
#include <type_traits>
|
||||||
|
|
||||||
|
// The Serializer class encodes data into a memory buffer so it can be later
|
||||||
|
// reconstructed. This is mainly used to send objects back-and-forth between
|
||||||
|
// the CPU and GPU.
|
||||||
|
|
||||||
|
namespace Eigen {
|
||||||
|
|
||||||
|
/**
|
||||||
|
* Serializes an object to a memory buffer.
|
||||||
|
*
|
||||||
|
* Useful for transfering data (e.g. back-and-forth to a device).
|
||||||
|
*/
|
||||||
|
template<typename T, typename EnableIf = void>
|
||||||
|
class Serializer;
|
||||||
|
|
||||||
|
// Specialization for POD types.
|
||||||
|
template<typename T>
|
||||||
|
class Serializer<T, typename std::enable_if<
|
||||||
|
std::is_trivial<T>::value
|
||||||
|
&& std::is_standard_layout<T>::value>::type > {
|
||||||
|
public:
|
||||||
|
|
||||||
|
/**
|
||||||
|
* Determines the required size of the serialization buffer for a value.
|
||||||
|
*
|
||||||
|
* \param value the value to serialize.
|
||||||
|
* \return the required size.
|
||||||
|
*/
|
||||||
|
EIGEN_DEVICE_FUNC size_t size(const T& value) const {
|
||||||
|
return sizeof(value);
|
||||||
|
}
|
||||||
|
|
||||||
|
/**
|
||||||
|
* Serializes a value to a byte buffer.
|
||||||
|
* \param dest the destination buffer.
|
||||||
|
* \param T the value to serialize.
|
||||||
|
* \return the next memory address past the end of the serialized data.
|
||||||
|
*/
|
||||||
|
EIGEN_DEVICE_FUNC uint8_t* serialize(uint8_t* dest, const T& value) {
|
||||||
|
EIGEN_USING_STD(memcpy)
|
||||||
|
memcpy(dest, &value, sizeof(value));
|
||||||
|
return dest + sizeof(value);
|
||||||
|
}
|
||||||
|
|
||||||
|
/**
|
||||||
|
* Deserializes a value from a byte buffer.
|
||||||
|
* \param src the source buffer.
|
||||||
|
* \param value the value to populate.
|
||||||
|
* \return the next unprocessed memory address.
|
||||||
|
*/
|
||||||
|
EIGEN_DEVICE_FUNC uint8_t* deserialize(uint8_t* src, T& value) const {
|
||||||
|
EIGEN_USING_STD(memcpy)
|
||||||
|
memcpy(&value, src, sizeof(value));
|
||||||
|
return src + sizeof(value);
|
||||||
|
}
|
||||||
|
};
|
||||||
|
|
||||||
|
// Specialization for DenseBase.
|
||||||
|
// Serializes [rows, cols, data...].
|
||||||
|
template<typename Derived>
|
||||||
|
class Serializer<DenseBase<Derived>, void> {
|
||||||
|
public:
|
||||||
|
typedef typename Derived::Scalar Scalar;
|
||||||
|
|
||||||
|
struct Header {
|
||||||
|
typename Derived::Index rows;
|
||||||
|
typename Derived::Index cols;
|
||||||
|
};
|
||||||
|
|
||||||
|
EIGEN_DEVICE_FUNC size_t size(const Derived& value) const {
|
||||||
|
return sizeof(Header) + sizeof(Scalar) * value.size();
|
||||||
|
}
|
||||||
|
|
||||||
|
EIGEN_DEVICE_FUNC uint8_t* serialize(uint8_t* dest, const Derived& value) {
|
||||||
|
const size_t header_bytes = sizeof(Header);
|
||||||
|
const size_t data_bytes = sizeof(Scalar) * value.size();
|
||||||
|
Header header = {value.rows(), value.cols()};
|
||||||
|
EIGEN_USING_STD(memcpy)
|
||||||
|
memcpy(dest, &header, header_bytes);
|
||||||
|
dest += header_bytes;
|
||||||
|
memcpy(dest, value.data(), data_bytes);
|
||||||
|
return dest + data_bytes;
|
||||||
|
}
|
||||||
|
|
||||||
|
EIGEN_DEVICE_FUNC uint8_t* deserialize(uint8_t* src, Derived& value) const {
|
||||||
|
const size_t header_bytes = sizeof(Header);
|
||||||
|
Header header;
|
||||||
|
EIGEN_USING_STD(memcpy)
|
||||||
|
memcpy(&header, src, header_bytes);
|
||||||
|
src += header_bytes;
|
||||||
|
value.resize(header.rows, header.cols);
|
||||||
|
const size_t data_bytes = sizeof(Scalar) * header.rows * header.cols;
|
||||||
|
memcpy(value.data(), src, data_bytes);
|
||||||
|
return src + data_bytes;
|
||||||
|
}
|
||||||
|
};
|
||||||
|
|
||||||
|
template<typename Scalar, int Rows, int Cols, int Options, int MaxRows, int MaxCols>
|
||||||
|
class Serializer<Matrix<Scalar, Rows, Cols, Options, MaxRows, MaxCols> > : public
|
||||||
|
Serializer<DenseBase<Matrix<Scalar, Rows, Cols, Options, MaxRows, MaxCols> > > {};
|
||||||
|
|
||||||
|
template<typename Scalar, int Rows, int Cols, int Options, int MaxRows, int MaxCols>
|
||||||
|
class Serializer<Array<Scalar, Rows, Cols, Options, MaxRows, MaxCols> > : public
|
||||||
|
Serializer<DenseBase<Array<Scalar, Rows, Cols, Options, MaxRows, MaxCols> > > {};
|
||||||
|
|
||||||
|
namespace internal {
|
||||||
|
|
||||||
|
// Recursive serialization implementation helper.
|
||||||
|
template<size_t N, typename... Types>
|
||||||
|
struct serialize_impl;
|
||||||
|
|
||||||
|
template<size_t N, typename T1, typename... Ts>
|
||||||
|
struct serialize_impl<N, T1, Ts...> {
|
||||||
|
using Serializer = Eigen::Serializer<typename std::decay<T1>::type>;
|
||||||
|
|
||||||
|
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
|
||||||
|
size_t serialize_size(const T1& value, const Ts&... args) {
|
||||||
|
Serializer serializer;
|
||||||
|
size_t size = serializer.size(value);
|
||||||
|
return size + serialize_impl<N-1, Ts...>::serialize_size(args...);
|
||||||
|
}
|
||||||
|
|
||||||
|
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
|
||||||
|
uint8_t* serialize(uint8_t* dest, const T1& value, const Ts&... args) {
|
||||||
|
Serializer serializer;
|
||||||
|
dest = serializer.serialize(dest, value);
|
||||||
|
return serialize_impl<N-1, Ts...>::serialize(dest, args...);
|
||||||
|
}
|
||||||
|
|
||||||
|
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
|
||||||
|
uint8_t* deserialize(uint8_t* src, T1& value, Ts&... args) {
|
||||||
|
Serializer serializer;
|
||||||
|
src = serializer.deserialize(src, value);
|
||||||
|
return serialize_impl<N-1, Ts...>::deserialize(src, args...);
|
||||||
|
}
|
||||||
|
};
|
||||||
|
|
||||||
|
// Base case.
|
||||||
|
template<>
|
||||||
|
struct serialize_impl<0> {
|
||||||
|
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
|
||||||
|
size_t serialize_size() { return 0; }
|
||||||
|
|
||||||
|
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
|
||||||
|
uint8_t* serialize(uint8_t* dest) { return dest; }
|
||||||
|
|
||||||
|
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
|
||||||
|
uint8_t* deserialize(uint8_t* src) { return src; }
|
||||||
|
};
|
||||||
|
|
||||||
|
} // namespace internal
|
||||||
|
|
||||||
|
|
||||||
|
/**
|
||||||
|
* Determine the buffer size required to serialize a set of values.
|
||||||
|
*
|
||||||
|
* \param args ... arguments to serialize in sequence.
|
||||||
|
* \return the total size of the required buffer.
|
||||||
|
*/
|
||||||
|
template<typename... Args>
|
||||||
|
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
|
||||||
|
size_t serialize_size(const Args&... args) {
|
||||||
|
return internal::serialize_impl<sizeof...(args), Args...>::serialize_size(args...);
|
||||||
|
}
|
||||||
|
|
||||||
|
/**
|
||||||
|
* Serialize a set of values to the byte buffer.
|
||||||
|
*
|
||||||
|
* \param dest output byte buffer.
|
||||||
|
* \param args ... arguments to serialize in sequence.
|
||||||
|
* \return the next address after all serialized values.
|
||||||
|
*/
|
||||||
|
template<typename... Args>
|
||||||
|
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
|
||||||
|
uint8_t* serialize(uint8_t* dest, const Args&... args) {
|
||||||
|
return internal::serialize_impl<sizeof...(args), Args...>::serialize(dest, args...);
|
||||||
|
}
|
||||||
|
|
||||||
|
/**
|
||||||
|
* Deserialize a set of values from the byte buffer.
|
||||||
|
*
|
||||||
|
* \param src input byte buffer.
|
||||||
|
* \param args ... arguments to deserialize in sequence.
|
||||||
|
* \return the next address after all parsed values.
|
||||||
|
*/
|
||||||
|
template<typename... Args>
|
||||||
|
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
|
||||||
|
uint8_t* deserialize(uint8_t* src, Args&... args) {
|
||||||
|
return internal::serialize_impl<sizeof...(args), Args...>::deserialize(src, args...);
|
||||||
|
}
|
||||||
|
|
||||||
|
} // namespace Eigen
|
||||||
|
|
||||||
|
#endif // EIGEN_SERIALIZER_H
|
@ -288,6 +288,7 @@ ei_add_test(blasutil)
|
|||||||
ei_add_test(random_matrix)
|
ei_add_test(random_matrix)
|
||||||
ei_add_test(initializer_list_construction)
|
ei_add_test(initializer_list_construction)
|
||||||
ei_add_test(diagonal_matrix_variadic_ctor)
|
ei_add_test(diagonal_matrix_variadic_ctor)
|
||||||
|
ei_add_test(serializer)
|
||||||
|
|
||||||
add_executable(bug1213 bug1213.cpp bug1213_main.cpp)
|
add_executable(bug1213 bug1213.cpp bug1213_main.cpp)
|
||||||
|
|
||||||
|
25
test/main.h
25
test/main.h
@ -391,6 +391,8 @@ inline void verify_impl(bool condition, const char *testname, const char *file,
|
|||||||
#define VERIFY_IS_NOT_MUCH_SMALLER_THAN(a, b) VERIFY(!test_isMuchSmallerThan(a, b))
|
#define VERIFY_IS_NOT_MUCH_SMALLER_THAN(a, b) VERIFY(!test_isMuchSmallerThan(a, b))
|
||||||
#define VERIFY_IS_APPROX_OR_LESS_THAN(a, b) VERIFY(test_isApproxOrLessThan(a, b))
|
#define VERIFY_IS_APPROX_OR_LESS_THAN(a, b) VERIFY(test_isApproxOrLessThan(a, b))
|
||||||
#define VERIFY_IS_NOT_APPROX_OR_LESS_THAN(a, b) VERIFY(!test_isApproxOrLessThan(a, b))
|
#define VERIFY_IS_NOT_APPROX_OR_LESS_THAN(a, b) VERIFY(!test_isApproxOrLessThan(a, b))
|
||||||
|
#define VERIFY_IS_CWISE_EQUAL(a, b) VERIFY(test_isCwiseApprox(a, b, true))
|
||||||
|
#define VERIFY_IS_CWISE_APPROX(a, b) VERIFY(test_isCwiseApprox(a, b, false))
|
||||||
|
|
||||||
#define VERIFY_IS_UNITARY(a) VERIFY(test_isUnitary(a))
|
#define VERIFY_IS_UNITARY(a) VERIFY(test_isUnitary(a))
|
||||||
|
|
||||||
@ -655,6 +657,29 @@ inline bool test_isUnitary(const MatrixBase<Derived>& m)
|
|||||||
return m.isUnitary(test_precision<typename internal::traits<Derived>::Scalar>());
|
return m.isUnitary(test_precision<typename internal::traits<Derived>::Scalar>());
|
||||||
}
|
}
|
||||||
|
|
||||||
|
// Checks component-wise, works with infs and nans.
|
||||||
|
template<typename Derived1, typename Derived2>
|
||||||
|
bool test_isCwiseApprox(const DenseBase<Derived1>& m1,
|
||||||
|
const DenseBase<Derived2>& m2,
|
||||||
|
bool exact) {
|
||||||
|
if (m1.rows() != m2.rows()) {
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
if (m1.cols() != m2.cols()) {
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
for (Index r = 0; r < m1.rows(); ++r) {
|
||||||
|
for (Index c = 0; c < m1.cols(); ++c) {
|
||||||
|
if (m1(r, c) != m2(r, c)
|
||||||
|
&& !((numext::isnan)(m1(r, c)) && (numext::isnan)(m2(r, c)))
|
||||||
|
&& (exact || !test_isApprox(m1(r, c), m2(r, c)))) {
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
template<typename T, typename U>
|
template<typename T, typename U>
|
||||||
bool test_is_equal(const T& actual, const U& expected, bool expect_equal)
|
bool test_is_equal(const T& actual, const U& expected, bool expect_equal)
|
||||||
{
|
{
|
||||||
|
108
test/serializer.cpp
Normal file
108
test/serializer.cpp
Normal file
@ -0,0 +1,108 @@
|
|||||||
|
// This file is part of Eigen, a lightweight C++ template library
|
||||||
|
// for linear algebra.
|
||||||
|
//
|
||||||
|
// Copyright (C) 2021 The Eigen Team
|
||||||
|
//
|
||||||
|
// This Source Code Form is subject to the terms of the Mozilla
|
||||||
|
// Public License v. 2.0. If a copy of the MPL was not distributed
|
||||||
|
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
||||||
|
|
||||||
|
#include "main.h"
|
||||||
|
|
||||||
|
#include <vector>
|
||||||
|
#include <Eigen/Core>
|
||||||
|
|
||||||
|
struct MyPodType {
|
||||||
|
double x;
|
||||||
|
int y;
|
||||||
|
float z;
|
||||||
|
};
|
||||||
|
|
||||||
|
// Plain-old-data serialization.
|
||||||
|
void test_pod_type() {
|
||||||
|
MyPodType initial = {1.3, 17, 1.9f};
|
||||||
|
MyPodType clone = {-1, -1, -1};
|
||||||
|
|
||||||
|
Eigen::Serializer<MyPodType> serializer;
|
||||||
|
|
||||||
|
// Determine required size.
|
||||||
|
size_t buffer_size = serializer.size(initial);
|
||||||
|
VERIFY_IS_EQUAL(buffer_size, sizeof(MyPodType));
|
||||||
|
|
||||||
|
// Serialize.
|
||||||
|
std::vector<uint8_t> buffer(buffer_size);
|
||||||
|
uint8_t* dest = serializer.serialize(buffer.data(), initial);
|
||||||
|
VERIFY_IS_EQUAL(dest - buffer.data(), buffer_size);
|
||||||
|
|
||||||
|
// Deserialize.
|
||||||
|
uint8_t* src = serializer.deserialize(buffer.data(), clone);
|
||||||
|
VERIFY_IS_EQUAL(src - buffer.data(), buffer_size);
|
||||||
|
VERIFY_IS_EQUAL(clone.x, initial.x);
|
||||||
|
VERIFY_IS_EQUAL(clone.y, initial.y);
|
||||||
|
VERIFY_IS_EQUAL(clone.z, initial.z);
|
||||||
|
}
|
||||||
|
|
||||||
|
// Matrix, Vector, Array
|
||||||
|
template<typename T>
|
||||||
|
void test_eigen_type(const T& type) {
|
||||||
|
const Index rows = type.rows();
|
||||||
|
const Index cols = type.cols();
|
||||||
|
|
||||||
|
const T initial = T::Random(rows, cols);
|
||||||
|
|
||||||
|
// Serialize.
|
||||||
|
Eigen::Serializer<T> serializer;
|
||||||
|
size_t buffer_size = serializer.size(initial);
|
||||||
|
std::vector<uint8_t> buffer(buffer_size);
|
||||||
|
uint8_t* dest = serializer.serialize(buffer.data(), initial);
|
||||||
|
VERIFY_IS_EQUAL(dest - buffer.data(), buffer_size);
|
||||||
|
|
||||||
|
// Deserialize.
|
||||||
|
T clone;
|
||||||
|
uint8_t* src = serializer.deserialize(buffer.data(), clone);
|
||||||
|
VERIFY_IS_EQUAL(src - buffer.data(), buffer_size);
|
||||||
|
VERIFY_IS_CWISE_EQUAL(clone, initial);
|
||||||
|
}
|
||||||
|
|
||||||
|
// Test a collection of dense types.
|
||||||
|
template<typename T1, typename T2, typename T3>
|
||||||
|
void test_dense_types(const T1& type1, const T2& type2, const T3& type3) {
|
||||||
|
|
||||||
|
// Make random inputs.
|
||||||
|
const T1 x1 = T1::Random(type1.rows(), type1.cols());
|
||||||
|
const T2 x2 = T2::Random(type2.rows(), type2.cols());
|
||||||
|
const T3 x3 = T3::Random(type3.rows(), type3.cols());
|
||||||
|
|
||||||
|
// Allocate buffer and serialize.
|
||||||
|
size_t buffer_size = Eigen::serialize_size(x1, x2, x3);
|
||||||
|
std::vector<uint8_t> buffer(buffer_size);
|
||||||
|
Eigen::serialize(buffer.data(), x1, x2, x3);
|
||||||
|
|
||||||
|
// Clone everything.
|
||||||
|
T1 y1;
|
||||||
|
T2 y2;
|
||||||
|
T3 y3;
|
||||||
|
Eigen::deserialize(buffer.data(), y1, y2, y3);
|
||||||
|
|
||||||
|
// Verify they equal.
|
||||||
|
VERIFY_IS_CWISE_EQUAL(y1, x1);
|
||||||
|
VERIFY_IS_CWISE_EQUAL(y2, x2);
|
||||||
|
VERIFY_IS_CWISE_EQUAL(y3, x3);
|
||||||
|
}
|
||||||
|
|
||||||
|
EIGEN_DECLARE_TEST(serializer)
|
||||||
|
{
|
||||||
|
CALL_SUBTEST( test_pod_type() );
|
||||||
|
|
||||||
|
for(int i = 0; i < g_repeat; i++) {
|
||||||
|
CALL_SUBTEST( test_eigen_type(Eigen::Array33f()) );
|
||||||
|
CALL_SUBTEST( test_eigen_type(Eigen::ArrayXd(10)) );
|
||||||
|
CALL_SUBTEST( test_eigen_type(Eigen::Vector3f()) );
|
||||||
|
CALL_SUBTEST( test_eigen_type(Eigen::Matrix4d()) );
|
||||||
|
CALL_SUBTEST( test_eigen_type(Eigen::MatrixXd(15, 17)) );
|
||||||
|
|
||||||
|
CALL_SUBTEST( test_dense_types( Eigen::Array33f(),
|
||||||
|
Eigen::ArrayXd(10),
|
||||||
|
Eigen::MatrixXd(15, 17)) );
|
||||||
|
}
|
||||||
|
}
|
Loading…
x
Reference in New Issue
Block a user