123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406 |
- #include <utility>
- #include "ceres/ceres.h"
- #include "glog/logging.h"
- using EigenDouble = Eigen::NumTraits<double>;
- using Mat = Eigen::MatrixXd;
- using Vec = Eigen::VectorXd;
- using Mat3 = Eigen::Matrix<double, 3, 3>;
- using Vec2 = Eigen::Matrix<double, 2, 1>;
- using MatX8 = Eigen::Matrix<double, Eigen::Dynamic, 8>;
- using Vec3 = Eigen::Vector3d;
- namespace {
- struct EstimateHomographyOptions {
-
-
- EstimateHomographyOptions() = default;
-
- int max_num_iterations{50};
-
-
-
-
-
-
-
-
- double expected_average_symmetric_distance{1e-16};
- };
- template <typename T>
- void SymmetricGeometricDistanceTerms(const Eigen::Matrix<T, 3, 3>& H,
- const Eigen::Matrix<T, 2, 1>& x1,
- const Eigen::Matrix<T, 2, 1>& x2,
- T forward_error[2],
- T backward_error[2]) {
- using Vec3 = Eigen::Matrix<T, 3, 1>;
- Vec3 x(x1(0), x1(1), T(1.0));
- Vec3 y(x2(0), x2(1), T(1.0));
- Vec3 H_x = H * x;
- Vec3 Hinv_y = H.inverse() * y;
- H_x /= H_x(2);
- Hinv_y /= Hinv_y(2);
- forward_error[0] = H_x(0) - y(0);
- forward_error[1] = H_x(1) - y(1);
- backward_error[0] = Hinv_y(0) - x(0);
- backward_error[1] = Hinv_y(1) - x(1);
- }
- double SymmetricGeometricDistance(const Mat3& H,
- const Vec2& x1,
- const Vec2& x2) {
- Vec2 forward_error, backward_error;
- SymmetricGeometricDistanceTerms<double>(
- H, x1, x2, forward_error.data(), backward_error.data());
- return forward_error.squaredNorm() + backward_error.squaredNorm();
- }
- template <typename T = double>
- class Homography2DNormalizedParameterization {
- public:
- using Parameters = Eigen::Matrix<T, 8, 1>;
- using Parameterized = Eigen::Matrix<T, 3, 3>;
-
- static void To(const Parameters& p, Parameterized* h) {
-
- *h << p(0), p(1), p(2),
- p(3), p(4), p(5),
- p(6), p(7), 1.0;
-
- }
-
- static void From(const Parameterized& h, Parameters* p) {
-
- *p << h(0, 0), h(0, 1), h(0, 2),
- h(1, 0), h(1, 1), h(1, 2),
- h(2, 0), h(2, 1);
-
- }
- };
- bool Homography2DFromCorrespondencesLinearEuc(const Mat& x1,
- const Mat& x2,
- Mat3* H,
- double expected_precision) {
- assert(2 == x1.rows());
- assert(4 <= x1.cols());
- assert(x1.rows() == x2.rows());
- assert(x1.cols() == x2.cols());
- const int64_t n = x1.cols();
- MatX8 L = Mat::Zero(n * 3, 8);
- Mat b = Mat::Zero(n * 3, 1);
- for (int64_t i = 0; i < n; ++i) {
- int64_t j = 3 * i;
- L(j, 0) = x1(0, i);
- L(j, 1) = x1(1, i);
- L(j, 2) = 1.0;
- L(j, 6) = -x2(0, i) * x1(0, i);
- L(j, 7) = -x2(0, i) * x1(1, i);
- b(j, 0) = x2(0, i);
- ++j;
- L(j, 3) = x1(0, i);
- L(j, 4) = x1(1, i);
- L(j, 5) = 1.0;
- L(j, 6) = -x2(1, i) * x1(0, i);
- L(j, 7) = -x2(1, i) * x1(1, i);
- b(j, 0) = x2(1, i);
-
-
- ++j;
- L(j, 0) = x2(1, i) * x1(0, i);
- L(j, 1) = x2(1, i) * x1(1, i);
- L(j, 2) = x2(1, i);
- L(j, 3) = -x2(0, i) * x1(0, i);
- L(j, 4) = -x2(0, i) * x1(1, i);
- L(j, 5) = -x2(0, i);
- }
-
- const Vec h = L.fullPivLu().solve(b);
- Homography2DNormalizedParameterization<double>::To(h, H);
- return (L * h).isApprox(b, expected_precision);
- }
- class HomographySymmetricGeometricCostFunctor {
- public:
- HomographySymmetricGeometricCostFunctor(Vec2 x, Vec2 y)
- : x_(std::move(x)), y_(std::move(y)) {}
- template <typename T>
- bool operator()(const T* homography_parameters, T* residuals) const {
- using Mat3 = Eigen::Matrix<T, 3, 3>;
- using Vec2 = Eigen::Matrix<T, 2, 1>;
- Mat3 H(homography_parameters);
- Vec2 x(T(x_(0)), T(x_(1)));
- Vec2 y(T(y_(0)), T(y_(1)));
- SymmetricGeometricDistanceTerms<T>(H, x, y, &residuals[0], &residuals[2]);
- return true;
- }
- const Vec2 x_;
- const Vec2 y_;
- };
- class TerminationCheckingCallback : public ceres::IterationCallback {
- public:
- TerminationCheckingCallback(const Mat& x1,
- const Mat& x2,
- const EstimateHomographyOptions& options,
- Mat3* H)
- : options_(options), x1_(x1), x2_(x2), H_(H) {}
- ceres::CallbackReturnType operator()(
- const ceres::IterationSummary& summary) override {
-
- if (!summary.step_is_successful) {
- return ceres::SOLVER_CONTINUE;
- }
-
- double average_distance = 0.0;
- for (int i = 0; i < x1_.cols(); i++) {
- average_distance +=
- SymmetricGeometricDistance(*H_, x1_.col(i), x2_.col(i));
- }
- average_distance /= x1_.cols();
- if (average_distance <= options_.expected_average_symmetric_distance) {
- return ceres::SOLVER_TERMINATE_SUCCESSFULLY;
- }
- return ceres::SOLVER_CONTINUE;
- }
- private:
- const EstimateHomographyOptions& options_;
- const Mat& x1_;
- const Mat& x2_;
- Mat3* H_;
- };
- bool EstimateHomography2DFromCorrespondences(
- const Mat& x1,
- const Mat& x2,
- const EstimateHomographyOptions& options,
- Mat3* H) {
- assert(2 == x1.rows());
- assert(4 <= x1.cols());
- assert(x1.rows() == x2.rows());
- assert(x1.cols() == x2.cols());
-
-
- Homography2DFromCorrespondencesLinearEuc(
- x1, x2, H, EigenDouble::dummy_precision());
- LOG(INFO) << "Estimated matrix after algebraic estimation:\n" << *H;
-
- ceres::Problem problem;
- for (int i = 0; i < x1.cols(); i++) {
- auto* homography_symmetric_geometric_cost_function =
- new HomographySymmetricGeometricCostFunctor(x1.col(i), x2.col(i));
- problem.AddResidualBlock(
- new ceres::AutoDiffCostFunction<HomographySymmetricGeometricCostFunctor,
- 4,
- 9>(
- homography_symmetric_geometric_cost_function),
- nullptr,
- H->data());
- }
-
- ceres::Solver::Options solver_options;
- solver_options.linear_solver_type = ceres::DENSE_QR;
- solver_options.max_num_iterations = options.max_num_iterations;
- solver_options.update_state_every_iteration = true;
-
- TerminationCheckingCallback callback(x1, x2, options, H);
- solver_options.callbacks.push_back(&callback);
-
- ceres::Solver::Summary summary;
- ceres::Solve(solver_options, &problem, &summary);
- LOG(INFO) << "Summary:\n" << summary.FullReport();
- LOG(INFO) << "Final refined matrix:\n" << *H;
- return summary.IsSolutionUsable();
- }
- }
- int main(int argc, char** argv) {
- google::InitGoogleLogging(argv[0]);
- Mat x1(2, 100);
- for (int i = 0; i < x1.cols(); ++i) {
- x1(0, i) = rand() % 1024;
- x1(1, i) = rand() % 1024;
- }
- Mat3 homography_matrix;
-
-
- homography_matrix << 1.243715, -0.461057, -111.964454,
- 0.0, 0.617589, -192.379252,
- 0.0, -0.000983, 1.0;
-
- Mat x2 = x1;
- for (int i = 0; i < x2.cols(); ++i) {
- Vec3 homogeneous_x1 = Vec3(x1(0, i), x1(1, i), 1.0);
- Vec3 homogeneous_x2 = homography_matrix * homogeneous_x1;
- x2(0, i) = homogeneous_x2(0) / homogeneous_x2(2);
- x2(1, i) = homogeneous_x2(1) / homogeneous_x2(2);
-
- x2(0, i) += static_cast<double>(rand() % 1000) / 5000.0;
- x2(1, i) += static_cast<double>(rand() % 1000) / 5000.0;
- }
- Mat3 estimated_matrix;
- EstimateHomographyOptions options;
- options.expected_average_symmetric_distance = 0.02;
- EstimateHomography2DFromCorrespondences(x1, x2, options, &estimated_matrix);
-
- estimated_matrix /= estimated_matrix(2, 2);
- std::cout << "Original matrix:\n" << homography_matrix << "\n";
- std::cout << "Estimated matrix:\n" << estimated_matrix << "\n";
- return EXIT_SUCCESS;
- }
|