Blame examples/simple_bundle_adjuster.cc

Packit ea1746
// Ceres Solver - A fast non-linear least squares minimizer
Packit ea1746
// Copyright 2015 Google Inc. All rights reserved.
Packit ea1746
// http://ceres-solver.org/
Packit ea1746
//
Packit ea1746
// Redistribution and use in source and binary forms, with or without
Packit ea1746
// modification, are permitted provided that the following conditions are met:
Packit ea1746
//
Packit ea1746
// * Redistributions of source code must retain the above copyright notice,
Packit ea1746
//   this list of conditions and the following disclaimer.
Packit ea1746
// * Redistributions in binary form must reproduce the above copyright notice,
Packit ea1746
//   this list of conditions and the following disclaimer in the documentation
Packit ea1746
//   and/or other materials provided with the distribution.
Packit ea1746
// * Neither the name of Google Inc. nor the names of its contributors may be
Packit ea1746
//   used to endorse or promote products derived from this software without
Packit ea1746
//   specific prior written permission.
Packit ea1746
//
Packit ea1746
// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
Packit ea1746
// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
Packit ea1746
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
Packit ea1746
// ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
Packit ea1746
// LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
Packit ea1746
// CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
Packit ea1746
// SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
Packit ea1746
// INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
Packit ea1746
// CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
Packit ea1746
// ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
Packit ea1746
// POSSIBILITY OF SUCH DAMAGE.
Packit ea1746
//
Packit ea1746
// Author: keir@google.com (Keir Mierle)
Packit ea1746
//
Packit ea1746
// A minimal, self-contained bundle adjuster using Ceres, that reads
Packit ea1746
// files from University of Washington' Bundle Adjustment in the Large dataset:
Packit ea1746
// http://grail.cs.washington.edu/projects/bal
Packit ea1746
//
Packit ea1746
// This does not use the best configuration for solving; see the more involved
Packit ea1746
// bundle_adjuster.cc file for details.
Packit ea1746
Packit ea1746
#include <cmath>
Packit ea1746
#include <cstdio>
Packit ea1746
#include <iostream>
Packit ea1746
Packit ea1746
#include "ceres/ceres.h"
Packit ea1746
#include "ceres/rotation.h"
Packit ea1746
Packit ea1746
// Read a Bundle Adjustment in the Large dataset.
Packit ea1746
class BALProblem {
Packit ea1746
 public:
Packit ea1746
  ~BALProblem() {
Packit ea1746
    delete[] point_index_;
Packit ea1746
    delete[] camera_index_;
Packit ea1746
    delete[] observations_;
Packit ea1746
    delete[] parameters_;
Packit ea1746
  }
Packit ea1746
Packit ea1746
  int num_observations()       const { return num_observations_;               }
Packit ea1746
  const double* observations() const { return observations_;                   }
Packit ea1746
  double* mutable_cameras()          { return parameters_;                     }
Packit ea1746
  double* mutable_points()           { return parameters_  + 9 * num_cameras_; }
Packit ea1746
Packit ea1746
  double* mutable_camera_for_observation(int i) {
Packit ea1746
    return mutable_cameras() + camera_index_[i] * 9;
Packit ea1746
  }
Packit ea1746
  double* mutable_point_for_observation(int i) {
Packit ea1746
    return mutable_points() + point_index_[i] * 3;
Packit ea1746
  }
Packit ea1746
Packit ea1746
  bool LoadFile(const char* filename) {
Packit ea1746
    FILE* fptr = fopen(filename, "r");
Packit ea1746
    if (fptr == NULL) {
Packit ea1746
      return false;
Packit ea1746
    };
Packit ea1746
Packit ea1746
    FscanfOrDie(fptr, "%d", &num_cameras_);
Packit ea1746
    FscanfOrDie(fptr, "%d", &num_points_);
Packit ea1746
    FscanfOrDie(fptr, "%d", &num_observations_);
Packit ea1746
Packit ea1746
    point_index_ = new int[num_observations_];
Packit ea1746
    camera_index_ = new int[num_observations_];
Packit ea1746
    observations_ = new double[2 * num_observations_];
Packit ea1746
Packit ea1746
    num_parameters_ = 9 * num_cameras_ + 3 * num_points_;
Packit ea1746
    parameters_ = new double[num_parameters_];
Packit ea1746
Packit ea1746
    for (int i = 0; i < num_observations_; ++i) {
Packit ea1746
      FscanfOrDie(fptr, "%d", camera_index_ + i);
Packit ea1746
      FscanfOrDie(fptr, "%d", point_index_ + i);
Packit ea1746
      for (int j = 0; j < 2; ++j) {
Packit ea1746
        FscanfOrDie(fptr, "%lf", observations_ + 2*i + j);
Packit ea1746
      }
Packit ea1746
    }
Packit ea1746
Packit ea1746
    for (int i = 0; i < num_parameters_; ++i) {
Packit ea1746
      FscanfOrDie(fptr, "%lf", parameters_ + i);
Packit ea1746
    }
Packit ea1746
    return true;
Packit ea1746
  }
Packit ea1746
Packit ea1746
 private:
Packit ea1746
  template<typename T>
Packit ea1746
  void FscanfOrDie(FILE *fptr, const char *format, T *value) {
Packit ea1746
    int num_scanned = fscanf(fptr, format, value);
Packit ea1746
    if (num_scanned != 1) {
Packit ea1746
      LOG(FATAL) << "Invalid UW data file.";
Packit ea1746
    }
Packit ea1746
  }
Packit ea1746
Packit ea1746
  int num_cameras_;
Packit ea1746
  int num_points_;
Packit ea1746
  int num_observations_;
Packit ea1746
  int num_parameters_;
Packit ea1746
Packit ea1746
  int* point_index_;
Packit ea1746
  int* camera_index_;
Packit ea1746
  double* observations_;
Packit ea1746
  double* parameters_;
Packit ea1746
};
Packit ea1746
Packit ea1746
// Templated pinhole camera model for used with Ceres.  The camera is
Packit ea1746
// parameterized using 9 parameters: 3 for rotation, 3 for translation, 1 for
Packit ea1746
// focal length and 2 for radial distortion. The principal point is not modeled
Packit ea1746
// (i.e. it is assumed be located at the image center).
Packit ea1746
struct SnavelyReprojectionError {
Packit ea1746
  SnavelyReprojectionError(double observed_x, double observed_y)
Packit ea1746
      : observed_x(observed_x), observed_y(observed_y) {}
Packit ea1746
Packit ea1746
  template <typename T>
Packit ea1746
  bool operator()(const T* const camera,
Packit ea1746
                  const T* const point,
Packit ea1746
                  T* residuals) const {
Packit ea1746
    // camera[0,1,2] are the angle-axis rotation.
Packit ea1746
    T p[3];
Packit ea1746
    ceres::AngleAxisRotatePoint(camera, point, p);
Packit ea1746
Packit ea1746
    // camera[3,4,5] are the translation.
Packit ea1746
    p[0] += camera[3];
Packit ea1746
    p[1] += camera[4];
Packit ea1746
    p[2] += camera[5];
Packit ea1746
Packit ea1746
    // Compute the center of distortion. The sign change comes from
Packit ea1746
    // the camera model that Noah Snavely's Bundler assumes, whereby
Packit ea1746
    // the camera coordinate system has a negative z axis.
Packit ea1746
    T xp = - p[0] / p[2];
Packit ea1746
    T yp = - p[1] / p[2];
Packit ea1746
Packit ea1746
    // Apply second and fourth order radial distortion.
Packit ea1746
    const T& l1 = camera[7];
Packit ea1746
    const T& l2 = camera[8];
Packit ea1746
    T r2 = xp*xp + yp*yp;
Packit ea1746
    T distortion = 1.0 + r2  * (l1 + l2  * r2);
Packit ea1746
Packit ea1746
    // Compute final projected point position.
Packit ea1746
    const T& focal = camera[6];
Packit ea1746
    T predicted_x = focal * distortion * xp;
Packit ea1746
    T predicted_y = focal * distortion * yp;
Packit ea1746
Packit ea1746
    // The error is the difference between the predicted and observed position.
Packit ea1746
    residuals[0] = predicted_x - observed_x;
Packit ea1746
    residuals[1] = predicted_y - observed_y;
Packit ea1746
Packit ea1746
    return true;
Packit ea1746
  }
Packit ea1746
Packit ea1746
  // Factory to hide the construction of the CostFunction object from
Packit ea1746
  // the client code.
Packit ea1746
  static ceres::CostFunction* Create(const double observed_x,
Packit ea1746
                                     const double observed_y) {
Packit ea1746
    return (new ceres::AutoDiffCostFunction<SnavelyReprojectionError, 2, 9, 3>(
Packit ea1746
                new SnavelyReprojectionError(observed_x, observed_y)));
Packit ea1746
  }
Packit ea1746
Packit ea1746
  double observed_x;
Packit ea1746
  double observed_y;
Packit ea1746
};
Packit ea1746
Packit ea1746
int main(int argc, char** argv) {
Packit ea1746
  google::InitGoogleLogging(argv[0]);
Packit ea1746
  if (argc != 2) {
Packit ea1746
    std::cerr << "usage: simple_bundle_adjuster <bal_problem>\n";
Packit ea1746
    return 1;
Packit ea1746
  }
Packit ea1746
Packit ea1746
  BALProblem bal_problem;
Packit ea1746
  if (!bal_problem.LoadFile(argv[1])) {
Packit ea1746
    std::cerr << "ERROR: unable to open file " << argv[1] << "\n";
Packit ea1746
    return 1;
Packit ea1746
  }
Packit ea1746
Packit ea1746
  const double* observations = bal_problem.observations();
Packit ea1746
Packit ea1746
  // Create residuals for each observation in the bundle adjustment problem. The
Packit ea1746
  // parameters for cameras and points are added automatically.
Packit ea1746
  ceres::Problem problem;
Packit ea1746
  for (int i = 0; i < bal_problem.num_observations(); ++i) {
Packit ea1746
    // Each Residual block takes a point and a camera as input and outputs a 2
Packit ea1746
    // dimensional residual. Internally, the cost function stores the observed
Packit ea1746
    // image location and compares the reprojection against the observation.
Packit ea1746
Packit ea1746
    ceres::CostFunction* cost_function =
Packit ea1746
        SnavelyReprojectionError::Create(observations[2 * i + 0],
Packit ea1746
                                         observations[2 * i + 1]);
Packit ea1746
    problem.AddResidualBlock(cost_function,
Packit ea1746
                             NULL /* squared loss */,
Packit ea1746
                             bal_problem.mutable_camera_for_observation(i),
Packit ea1746
                             bal_problem.mutable_point_for_observation(i));
Packit ea1746
  }
Packit ea1746
Packit ea1746
  // Make Ceres automatically detect the bundle structure. Note that the
Packit ea1746
  // standard solver, SPARSE_NORMAL_CHOLESKY, also works fine but it is slower
Packit ea1746
  // for standard bundle adjustment problems.
Packit ea1746
  ceres::Solver::Options options;
Packit ea1746
  options.linear_solver_type = ceres::DENSE_SCHUR;
Packit ea1746
  options.minimizer_progress_to_stdout = true;
Packit ea1746
Packit ea1746
  ceres::Solver::Summary summary;
Packit ea1746
  ceres::Solve(options, &problem, &summary);
Packit ea1746
  std::cout << summary.FullReport() << "\n";
Packit ea1746
  return 0;
Packit ea1746
}