13 #include <Eigen/Geometry>
14 #include <Eigen/Sparse>
22 const Eigen::MatrixXd& A,
23 const Eigen::VectorXd& b,
25 bool check_symmetric ,
29 check_symmetric = check_symmetric || check_psd;
30 if (check_symmetric && !A.isApprox(A.transpose())) {
31 LogWarning(
"check_symmetric failed, empty vector will be returned");
32 return std::make_tuple(
false, Eigen::VectorXd::Zero(b.rows()));
36 double det = A.determinant();
37 if (
fabs(det) < 1e-6 || std::isnan(det) || std::isinf(det)) {
38 LogWarning(
"check_det failed, empty vector will be returned");
39 return std::make_tuple(
false, Eigen::VectorXd::Zero(b.rows()));
45 Eigen::LLT<Eigen::MatrixXd> A_llt(A);
46 if (A_llt.info() == Eigen::NumericalIssue) {
47 LogWarning(
"check_psd failed, empty vector will be returned");
48 return std::make_tuple(
false, Eigen::VectorXd::Zero(b.rows()));
52 Eigen::VectorXd x(b.size());
55 Eigen::SparseMatrix<double> A_sparse = A.sparseView();
57 Eigen::SimplicialCholesky<Eigen::SparseMatrix<double>> A_chol;
58 A_chol.compute(A_sparse);
59 if (A_chol.info() == Eigen::Success) {
61 if (A_chol.info() == Eigen::Success) {
63 return std::make_tuple(
true, std::move(x));
65 LogWarning(
"Cholesky solve failed, switched to dense solver");
68 LogWarning(
"Cholesky decompose failed, switched to dense solver");
72 x = A.ldlt().solve(b);
73 return std::make_tuple(
true, std::move(x));
77 Eigen::Matrix4d output;
79 output.block<3, 3>(0, 0) =
80 (Eigen::AngleAxisd(input(2), Eigen::Vector3d::UnitZ()) *
81 Eigen::AngleAxisd(input(1), Eigen::Vector3d::UnitY()) *
82 Eigen::AngleAxisd(input(0), Eigen::Vector3d::UnitX()))
84 output.block<3, 1>(0, 3) = input.block<3, 1>(3, 0);
90 Eigen::Matrix3d R = input.block<3, 3>(0, 0);
91 double sy = sqrt(R(0, 0) * R(0, 0) + R(1, 0) * R(1, 0));
93 output(0) = atan2(R(2, 1), R(2, 2));
94 output(1) = atan2(-R(2, 0), sy);
95 output(2) = atan2(R(1, 0), R(0, 0));
97 output(0) = atan2(-R(1, 2), R(1, 1));
98 output(1) = atan2(-R(2, 0), sy);
101 output.block<3, 1>(3, 0) = input.block<3, 1>(0, 3);
111 if (solution_exist) {
113 return std::make_tuple(solution_exist, std::move(extrinsic));
115 return std::make_tuple(
false, Eigen::Matrix4d::Identity());
118 std::tuple<bool, std::vector<Eigen::Matrix4d, Matrix4d_allocator>>
120 const Eigen::VectorXd& JTr) {
121 std::vector<Eigen::Matrix4d, Matrix4d_allocator> output_matrix_array;
122 output_matrix_array.clear();
123 if (JTJ.rows() != JTr.rows() || JTJ.cols() % 6 != 0) {
125 "[SolveJacobianSystemAndObtainExtrinsicMatrixArray] "
126 "Unsupported matrix format.");
127 return std::make_tuple(
false, std::move(output_matrix_array));
134 if (solution_exist) {
135 int nposes = (int)x.rows() / 6;
136 for (
int i = 0; i < nposes; i++) {
137 Eigen::Matrix4d extrinsic =
139 output_matrix_array.push_back(extrinsic);
141 return std::make_tuple(solution_exist, std::move(output_matrix_array));
143 return std::make_tuple(
false, std::move(output_matrix_array));
147 template <
typename MatType,
typename VecType>
149 std::function<
void(
int, VecType&,
double&,
double&)> f,
161 double r2_sum_private = 0.0;
162 JTJ_private.setZero();
163 JTr_private.setZero();
167 #pragma omp for nowait
168 for (
int i = 0; i < iteration_num; i++) {
170 JTJ_private.noalias() += J_r * w * J_r.transpose();
171 JTr_private.noalias() += J_r * w * r;
172 r2_sum_private += r * r;
174 #pragma omp critical(ComputeJTJandJTr)
178 r2_sum += r2_sum_private;
182 LogDebug(
"Residual : {:.2e} (# of elements : {:d})",
183 r2_sum / (
double)iteration_num, iteration_num);
185 return std::make_tuple(std::move(JTJ), std::move(JTr), r2_sum);
188 template <
typename MatType,
typename VecType>
192 std::vector<VecType, Eigen::aligned_allocator<VecType>>&,
193 std::vector<double>&,
194 std::vector<double>&)> f,
206 double r2_sum_private = 0.0;
207 JTJ_private.setZero();
208 JTr_private.setZero();
209 std::vector<double> r;
210 std::vector<double> w;
211 std::vector<VecType, Eigen::aligned_allocator<VecType>> J_r;
212 #pragma omp for nowait
213 for (
int i = 0; i < iteration_num; i++) {
215 for (
int j = 0; j < (int)r.size(); j++) {
216 JTJ_private.noalias() += J_r[j] * w[j] * J_r[j].transpose();
217 JTr_private.noalias() += J_r[j] * w[j] * r[j];
218 r2_sum_private += r[j] * r[j];
221 #pragma omp critical(ComputeJTJandJTr)
225 r2_sum += r2_sum_private;
229 LogDebug(
"Residual : {:.2e} (# of elements : {:d})",
230 r2_sum / (
double)iteration_num, iteration_num);
232 return std::make_tuple(std::move(JTJ), std::move(JTr), r2_sum);
238 int iteration_num,
bool verbose);
241 std::function<
void(
int,
242 std::vector<Eigen::Vector6d, Vector6d_allocator>&,
243 std::vector<double>&,
244 std::vector<double>&)> f,
245 int iteration_num,
bool verbose);
250 rot << 1, 0, 0, 0, std::cos(radians), -std::sin(radians), 0,
251 std::sin(radians), std::cos(radians);
257 rot << std::cos(radians), 0, std::sin(radians), 0, 1, 0, -std::sin(radians),
258 0, std::cos(radians);
264 rot << std::cos(radians), -std::sin(radians), 0, std::sin(radians),
265 std::cos(radians), 0, 0, 0, 1;
271 for (
int i = 0; i < 3; ++i) {
279 return Eigen::Vector3d(r, g, b) / 255.0;
286 template <
typename IdxType>
288 const std::vector<IdxType>& indices) {
289 if (indices.empty()) {
290 return Eigen::Matrix3d::Identity();
292 Eigen::Matrix3d covariance;
293 Eigen::Matrix<double, 9, 1> cumulants;
295 for (
const auto& idx : indices) {
297 cumulants(0) +=
point(0);
298 cumulants(1) +=
point(1);
299 cumulants(2) +=
point(2);
307 cumulants /= (double)indices.size();
308 covariance(0, 0) = cumulants(3) - cumulants(0) * cumulants(0);
309 covariance(1, 1) = cumulants(6) - cumulants(1) * cumulants(1);
310 covariance(2, 2) = cumulants(8) - cumulants(2) * cumulants(2);
311 covariance(0, 1) = cumulants(4) - cumulants(0) * cumulants(1);
312 covariance(1, 0) = covariance(0, 1);
313 covariance(0, 2) = cumulants(5) - cumulants(0) * cumulants(2);
314 covariance(2, 0) = covariance(0, 2);
315 covariance(1, 2) = cumulants(7) - cumulants(1) * cumulants(2);
316 covariance(2, 1) = covariance(1, 2);
320 template <
typename IdxType>
322 const std::vector<Eigen::Vector3d>&
points,
323 const std::vector<IdxType>& indices) {
324 Eigen::Vector3d mean;
325 Eigen::Matrix3d covariance;
326 Eigen::Matrix<double, 9, 1> cumulants;
328 for (
const auto& idx : indices) {
330 cumulants(0) +=
point(0);
331 cumulants(1) +=
point(1);
332 cumulants(2) +=
point(2);
340 cumulants /= (double)indices.size();
341 mean(0) = cumulants(0);
342 mean(1) = cumulants(1);
343 mean(2) = cumulants(2);
344 covariance(0, 0) = cumulants(3) - cumulants(0) * cumulants(0);
345 covariance(1, 1) = cumulants(6) - cumulants(1) * cumulants(1);
346 covariance(2, 2) = cumulants(8) - cumulants(2) * cumulants(2);
347 covariance(0, 1) = cumulants(4) - cumulants(0) * cumulants(1);
348 covariance(1, 0) = covariance(0, 1);
349 covariance(0, 2) = cumulants(5) - cumulants(0) * cumulants(2);
350 covariance(2, 0) = covariance(0, 2);
351 covariance(1, 2) = cumulants(7) - cumulants(1) * cumulants(2);
352 covariance(2, 1) = covariance(1, 2);
353 return std::make_tuple(mean, covariance);
357 const std::vector<Eigen::Vector3d>&
points,
358 const std::vector<size_t>& indices);
360 const std::vector<Eigen::Vector3d>&
points,
361 const std::vector<size_t>& indices);
363 const std::vector<Eigen::Vector3d>&
points,
364 const std::vector<int>& indices);
366 const std::vector<Eigen::Vector3d>&
points,
367 const std::vector<int>& indices);
369 template <
typename RealType,
typename IdxType>
371 const RealType*
const points,
const std::vector<IdxType>& indices) {
372 Eigen::Vector3d mean;
373 Eigen::Matrix3d covariance;
374 Eigen::Matrix<double, 9, 1> cumulants;
376 for (
const auto& idx : indices) {
377 const Eigen::Map<const Eigen::Matrix<RealType, 3, 1>>
point(
points +
379 cumulants(0) +=
point(0);
380 cumulants(1) +=
point(1);
381 cumulants(2) +=
point(2);
389 cumulants /= (double)indices.size();
390 mean(0) = cumulants(0);
391 mean(1) = cumulants(1);
392 mean(2) = cumulants(2);
393 covariance(0, 0) = cumulants(3) - cumulants(0) * cumulants(0);
394 covariance(1, 1) = cumulants(6) - cumulants(1) * cumulants(1);
395 covariance(2, 2) = cumulants(8) - cumulants(2) * cumulants(2);
396 covariance(0, 1) = cumulants(4) - cumulants(0) * cumulants(1);
397 covariance(1, 0) = covariance(0, 1);
398 covariance(0, 2) = cumulants(5) - cumulants(0) * cumulants(2);
399 covariance(2, 0) = covariance(0, 2);
400 covariance(1, 2) = cumulants(7) - cumulants(1) * cumulants(2);
401 covariance(2, 1) = covariance(1, 2);
402 return std::make_tuple(mean, covariance);
406 const float*
const points,
const std::vector<size_t>& indices);
408 const double*
const points,
const std::vector<size_t>& indices);
410 const float*
const points,
const std::vector<int>& indices);
412 const double*
const points,
const std::vector<int>& indices);
415 Eigen::Matrix3d skew;
417 skew << 0, -vec.z(), vec.y(),
418 vec.z(), 0, -vec.x(),
419 -vec.y(), vec.x(), 0;
__host__ __device__ float2 fabs(float2 v)
Eigen::Matrix< double, 6, 1 > Vector6d
Eigen::Matrix< double, 6, 6 > Matrix6d
Extending Eigen namespace by adding frequently used matrix type.
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)
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::Matrix3d SkewMatrix(const Eigen::Vector3d &vec)
Genretate a skew-symmetric matrix from a vector 3x1.
Eigen::Vector3uint8 ColorToUint8(const Eigen::Vector3d &color)
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)
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)
Generic file read and write utility for python interface.