8 #ifndef CLOUDVIEWER_EIGEN_H
9 #define CLOUDVIEWER_EIGEN_H
14 #define CPP_VERSION _MSVC_LANG
16 #define CPP_VERSION __cplusplus
23 #ifndef EIGEN_HAS_CXX17_OVERALIGN
25 #if CPP_VERSION >= 201703L
26 #define EIGEN_HAS_CXX17_OVERALIGN 1
32 #include <Eigen/Geometry>
33 #include <Eigen/StdVector>
38 #include <unordered_map>
42 #if !EIGEN_VERSION_AT_LEAST(3, 4, 0) || CPP_VERSION < 201703L
44 #include <initializer_list>
45 #ifndef EIGEN_ALIGNED_ALLOCATOR
46 #define EIGEN_ALIGNED_ALLOCATOR Eigen::aligned_allocator
52 #define EIGEN_DEFINE_STL_VECTOR_SPECIALIZATION_CUSTOM(...) \
55 class vector<__VA_ARGS__, std::allocator<__VA_ARGS__>> \
56 : public vector<__VA_ARGS__, EIGEN_ALIGNED_ALLOCATOR<__VA_ARGS__>> { \
57 typedef vector<__VA_ARGS__, EIGEN_ALIGNED_ALLOCATOR<__VA_ARGS__>> \
61 typedef __VA_ARGS__ value_type; \
62 typedef vector_base::allocator_type allocator_type; \
63 typedef vector_base::size_type size_type; \
64 typedef vector_base::iterator iterator; \
65 explicit vector(const allocator_type& a = allocator_type()) \
67 template <typename InputIterator> \
68 vector(InputIterator first, \
70 const allocator_type& a = allocator_type()) \
71 : vector_base(first, last, a) {} \
72 vector(const vector& c) : vector_base(c) {} \
73 explicit vector(size_type num, const value_type& val = value_type()) \
74 : vector_base(num, val) {} \
75 vector(iterator start, iterator end) : vector_base(start, end) {} \
76 vector& operator=(const vector& x) { \
77 vector_base::operator=(x); \
80 vector(initializer_list<__VA_ARGS__> list) \
81 : vector_base(list.begin(), list.end()) {} \
112 typedef Eigen::Matrix<double, 6, 6, Eigen::DontAlign>
Matrix6d_u;
113 typedef Eigen::Matrix<double, 4, 4, Eigen::DontAlign>
Matrix4d_u;
114 typedef Eigen::Matrix<double, 3, 1, Eigen::DontAlign>
Vector3d_u;
115 typedef Eigen::Matrix<float, 3, 1, Eigen::DontAlign>
Vector3f_u;
116 typedef Eigen::Matrix<double, 4, 1, Eigen::DontAlign>
Vector4d_u;
117 typedef Eigen::Matrix<float, 4, 1, Eigen::DontAlign>
Vector4f_u;
153 const Eigen::VectorXd& b,
154 bool prefer_sparse =
false,
155 bool check_symmetric =
false,
156 bool check_det =
false,
157 bool check_psd =
false);
169 std::tuple<bool, std::vector<Eigen::Matrix4d, Matrix4d_allocator>>
171 const Eigen::MatrixXd& JTJ,
const Eigen::VectorXd& JTr);
178 template <
typename MatType,
typename VecType>
182 bool verbose =
true);
189 template <
typename MatType,
typename VecType>
193 std::vector<VecType, Eigen::aligned_allocator<VecType>>&,
194 std::vector<double>&,
195 std::vector<double>&)> f,
197 bool verbose =
true);
211 template <
typename IdxType>
214 const std::vector<IdxType>& indices);
217 template <
typename IdxType>
220 const std::vector<IdxType>& indices);
228 template <
typename RealType,
typename IdxType>
231 const std::vector<IdxType>& indices);
#define EIGEN_DEFINE_STL_VECTOR_SPECIALIZATION_CUSTOM(...)
Eigen::Matrix< double, 6, 1 > Vector6d
Eigen::Matrix< double, 4, 1, Eigen::DontAlign > Vector4d_u
Eigen::Matrix< double, 6, 6, Eigen::DontAlign > Matrix6d_u
Eigen::Matrix< double, 4, 4, Eigen::DontAlign > Matrix4d_u
Eigen::Matrix< float, 4, 1, Eigen::DontAlign > Vector4f_u
Eigen::Matrix< float, 3, 1, Eigen::DontAlign > Vector3f_u
Eigen::Matrix< double, 6, 6 > Matrix6d
Extending Eigen namespace by adding frequently used matrix type.
Eigen::Matrix< double, 3, 1, Eigen::DontAlign > Vector3d_u
Eigen::Matrix< uint8_t, 3, 1 > Vector3uint8
std::tuple< Eigen::Vector3d, Eigen::Matrix3d > ComputeMeanAndCovariance(const std::vector< Eigen::Vector3d > &points, const std::vector< IdxType > &indices)
Function to compute the mean and covariance matrix of a set of points.
std::tuple< bool, std::vector< Eigen::Matrix4d, Matrix4d_allocator > > SolveJacobianSystemAndObtainExtrinsicMatrixArray(const Eigen::MatrixXd &JTJ, const Eigen::VectorXd &JTr)
Eigen::Matrix3d RotationMatrixX(double radians)
Eigen::aligned_allocator< Eigen::Vector4i > Vector4i_allocator
Eigen::aligned_allocator< Eigen::Vector6d > Vector6d_allocator
std::tuple< bool, Eigen::VectorXd > SolveLinearSystemPSD(const Eigen::MatrixXd &A, const Eigen::VectorXd &b, bool prefer_sparse=false, bool check_symmetric=false, bool check_det=false, bool check_psd=false)
Function to solve Ax=b.
Eigen::Matrix3d RotationMatrixY(double radians)
Eigen::Matrix3d ComputeCovariance(const std::vector< Eigen::Vector3d > &points, const std::vector< IdxType > &indices)
Function to compute the covariance matrix of a set of points.
Eigen::Vector6d TransformMatrix4dToVector6d(const Eigen::Matrix4d &input)
Eigen::aligned_allocator< Eigen::Matrix4d > Matrix4d_allocator
Eigen::Matrix3d SkewMatrix(const Eigen::Vector3d &vec)
Genretate a skew-symmetric matrix from a vector 3x1.
Eigen::Vector3uint8 ColorToUint8(const Eigen::Vector3d &color)
Eigen::aligned_allocator< Eigen::Vector4d > Vector4d_allocator
Eigen::Matrix4d TransformVector6dToMatrix4d(const Eigen::Vector6d &input)
Eigen::Vector3d ColorToDouble(uint8_t r, uint8_t g, uint8_t b)
Color conversion from uint8_t 0-255 to double [0,1].
Eigen::Matrix3d RotationMatrixZ(double radians)
Eigen::aligned_allocator< Eigen::Vector2d > Vector2d_allocator
std::tuple< bool, Eigen::Matrix4d > SolveJacobianSystemAndObtainExtrinsicMatrix(const Eigen::Matrix6d &JTJ, const Eigen::Vector6d &JTr)
std::tuple< MatType, VecType, double > ComputeJTJandJTr(std::function< void(int, VecType &, double &, double &)> f, int iteration_num, bool verbose=true)
Eigen::aligned_allocator< Eigen::Matrix6d > Matrix6d_allocator
Eigen::aligned_allocator< Eigen::Vector3uint8 > Vector3uint8_allocator
Generic file read and write utility for python interface.