Berkeley SfM
Public Member Functions | Private Member Functions | Private Attributes | List of all members
bsfm::PoseEstimator2D3D Class Reference

#include <pose_estimator_2d3d.h>

Public Member Functions

 PoseEstimator2D3D ()
 
 ~PoseEstimator2D3D ()
 
bool Initialize (const FeatureList &points_2d, const Point3DList &points_3d, const CameraIntrinsics &intrinsics)
 
bool Solve (Pose &camera_pose)
 

Private Member Functions

bool ComputeInitialSolution (Matrix34d &initial_solution)
 
bool OptimizeSolution (Matrix34d &solution)
 
bool ExtractPose (const Matrix34d &P, Pose &pose)
 

Private Attributes

FeatureList points_2d_
 
Point3DList points_3d_
 
CameraIntrinsics intrinsics_
 
Matrix3d T_
 
Matrix4d U_
 

Detailed Description

Definition at line 76 of file pose_estimator_2d3d.h.

Constructor & Destructor Documentation

bsfm::PoseEstimator2D3D::PoseEstimator2D3D ( )

Definition at line 49 of file pose_estimator_2d3d.cpp.

50  : T_(Matrix3d::Identity()), U_(Matrix4d::Identity()) {}
bsfm::PoseEstimator2D3D::~PoseEstimator2D3D ( )

Definition at line 52 of file pose_estimator_2d3d.cpp.

52 {}

Member Function Documentation

bool bsfm::PoseEstimator2D3D::ComputeInitialSolution ( Matrix34d initial_solution)
private

Definition at line 122 of file pose_estimator_2d3d.cpp.

122  {
123  // Use Eq. 7.2 from H&Z: Multiple-View Geometry to determine an
124  // over-determined solution for the camera projection matrix P.
125 
126  // First build the A matrix, which is 2n x 12.
127  MatrixXd A;
128  A.resize(points_2d_.size() * 2, 12);
129  for (size_t ii = 0; ii < points_2d_.size(); ii++) {
130  // Get (u, v) image-space point.
131  const double u = points_2d_[ii].u_;
132  const double v = points_2d_[ii].v_;
133 
134  // Get (x, y, z) world-space point.
135  const double x = points_3d_[ii].X();
136  const double y = points_3d_[ii].Y();
137  const double z = points_3d_[ii].Z();
138 
139  // First 4 columns, top row of block.
140  A(2*ii+0, 0) = 0;
141  A(2*ii+0, 1) = 0;
142  A(2*ii+0, 2) = 0;
143  A(2*ii+0, 3) = 0;
144 
145  // Second 4 columns, top row of block.
146  A(2*ii+0, 4) = -x;
147  A(2*ii+0, 5) = -y;
148  A(2*ii+0, 6) = -z;
149  A(2*ii+0, 7) = -1;
150 
151  // Third 4 columns, top row of block.
152  A(2*ii+0, 8) = v*x;
153  A(2*ii+0, 9) = v*y;
154  A(2*ii+0, 10) = v*z;
155  A(2*ii+0, 11) = v;
156 
157  // First 4 columns, bottom row of block.
158  A(2*ii+1, 0) = x;
159  A(2*ii+1, 1) = y;
160  A(2*ii+1, 2) = z;
161  A(2*ii+1, 3) = 1;
162 
163  // Second 4 columns, bottom row of block.
164  A(2*ii+1, 4) = 0;
165  A(2*ii+1, 5) = 0;
166  A(2*ii+1, 6) = 0;
167  A(2*ii+1, 7) = 0;
168 
169  // Third 4 columns, bottom row of block.
170  A(2*ii+1, 8) = -u*x;
171  A(2*ii+1, 9) = -u*y;
172  A(2*ii+1, 10) = -u*z;
173  A(2*ii+1, 11) = -u;
174  }
175 
176  // Get svd(A). Save some time and compute a thin U. We still need a full V.
177  Eigen::JacobiSVD<Eigen::MatrixXd> svd;
178  svd.compute(A, Eigen::ComputeThinU | Eigen::ComputeFullV);
179  if (!svd.computeV()) {
180  VLOG(1) << "Failed to compute a singular value decomposition of A matrix.";
181  return false;
182  }
183 
184  // Get the projection matrix elements from the SVD decomposition.
185  const VectorXd P_vec = svd.matrixV().col(11).normalized();
186 
187  // Reshape the vector into the initial solution.
188  initial_solution.row(0) = P_vec.topRows(4).transpose();
189  initial_solution.row(1) = P_vec.middleRows(4, 4).transpose();
190  initial_solution.row(2) = P_vec.bottomRows(4).transpose();
191 
192  return true;
193 }
bool bsfm::PoseEstimator2D3D::ExtractPose ( const Matrix34d P,
Pose pose 
)
private

Definition at line 222 of file pose_estimator_2d3d.cpp.

222  {
223  // Un-normalize the projection matrix.
224  const Matrix34d P_unnormalized = T_.inverse() * P * U_;
225 
226  // Extract camera extrinsics matrix.
227  Matrix34d Rt = intrinsics_.Kinv() * P_unnormalized;
228 
229  // [R|t] is only determined up to scale. To get this scale, note that det(R)
230  // must equal 1. Also note that for an nxn matrix, c^n*det(R) = det(cR).
231  // Use this property to scale our matrix.
232  double det = Rt.block(0, 0, 3, 3).determinant();
233  if (std::abs(det) < 1e-16) {
234  LOG(WARNING) << "Computed rotation has a determinant of 0.";
235  return false;
236  }
237 
238  // Make sure the determinant is positive.
239  if (det < 0.0) {
240  det *= -1.0;
241  Rt *= -1.0;
242  }
243 
244  // Normalize the rotation and translation.
245  Rt *= std::pow(1.0 / det, 1.0 / 3.0);
246 
247  // Initialize the output pose from the rotation and translation blocks.
248  pose = Pose(Rt);
249  return true;
250 }
::Eigen::Matrix< double, 3, 4 > Matrix34d
Definition: types.h:87
bool bsfm::PoseEstimator2D3D::Initialize ( const FeatureList points_2d,
const Point3DList points_3d,
const CameraIntrinsics intrinsics 
)

Definition at line 54 of file pose_estimator_2d3d.cpp.

56  {
57  if (points_2d.size() != points_3d.size()) {
58  VLOG(1) << "Inputs 'points_2d' and 'points_3d' do not contain the "
59  "same number of elements.";
60  return false;
61  }
62 
63  // Copy camera intrinsics.
64  intrinsics_ = intrinsics;
65 
66  // Normalize the 3D points.
67  T_ = ComputeNormalization(points_2d);
68  U_ = ComputeNormalization(points_3d);
69 
70  points_2d_.clear();
71  for (size_t ii = 0; ii < points_2d.size(); ++ii) {
72  // Get (u, v) image-space point and normalize.
73  double u = points_2d[ii].u_;
74  double v = points_2d[ii].v_;
75  u = T_(0, 0) * u + T_(0, 2);
76  v = T_(1, 1) * v + T_(1, 2);
77 
78  points_2d_.push_back(Feature(u, v));
79  }
80 
81  points_3d_.clear();
82  for (size_t ii = 0; ii < points_3d.size(); ++ii) {
83  // Get (x, y, z) world-space point and normalize.
84  double x = points_3d[ii].X();
85  double y = points_3d[ii].Y();
86  double z = points_3d[ii].Z();
87  x = U_(0, 0) * x + U_(0, 3);
88  y = U_(1, 1) * y + U_(1, 3);
89  z = U_(2, 2) * z + U_(2, 3);
90 
91  points_3d_.push_back(Point3D(x, y, z));
92  }
93 
94  return true;
95 }
Matrix3d ComputeNormalization(const FeatureMatchList &matched_features, bool use_feature_set1)
bool bsfm::PoseEstimator2D3D::OptimizeSolution ( Matrix34d solution)
private

Definition at line 195 of file pose_estimator_2d3d.cpp.

195  {
196  // Create the non-linear least squares problem and cost function.
197  ceres::Problem problem;
198 
199  // Add a term to the cost function for each 2D<-->3D correspondence.
200  Matrix34d optimized_solution = solution;
201  for (size_t ii = 0; ii < points_2d_.size(); ++ii) {
202  problem.AddResidualBlock(
204  NULL, /* squared loss */
205  optimized_solution.data());
206  }
207 
208  // Solve the non-linear least squares problem to get the projection matrix.
209  ceres::Solver::Summary summary;
210  ceres::Solver::Options options;
211  options.trust_region_strategy_type = ceres::LEVENBERG_MARQUARDT;
212  ceres::Solve(options, &problem, &summary);
213 
214  // Store the solved variable back in 'solution'.
215  if (summary.IsSolutionUsable()) {
216  solution = optimized_solution.normalized();
217  }
218 
219  return summary.IsSolutionUsable();
220 }
static ceres::CostFunction * Create(const Feature &x, const Point3D &X)
Definition: cost_functors.h:96
::Eigen::Matrix< double, 3, 4 > Matrix34d
Definition: types.h:87
bool bsfm::PoseEstimator2D3D::Solve ( Pose camera_pose)

Definition at line 97 of file pose_estimator_2d3d.cpp.

97  {
98  // Get an initial P.
99  Matrix34d P;
100  if (!ComputeInitialSolution(P)) {
101  VLOG(1) << "Failed to compute an initial solution for P.";
102  return false;
103  }
104 
105  // Refine P with non-linear optimization.
106  Matrix34d P_opt = P;
107 
108  if (!OptimizeSolution(P_opt)) {
109  VLOG(1) << "Failed to optimize P. Continuing using the initial solution.";
110  P_opt = P;
111  }
112 
113  // Get the camera pose from the computed projection matrix.
114  if (!ExtractPose(P_opt, camera_pose)) {
115  VLOG(1) << "Computed rotation is non-invertible.";
116  return false;
117  }
118 
119  return true;
120 }
bool ExtractPose(const Matrix34d &P, Pose &pose)
bool OptimizeSolution(Matrix34d &solution)
::Eigen::Matrix< double, 3, 4 > Matrix34d
Definition: types.h:87
bool ComputeInitialSolution(Matrix34d &initial_solution)

Member Data Documentation

CameraIntrinsics bsfm::PoseEstimator2D3D::intrinsics_
private

Definition at line 112 of file pose_estimator_2d3d.h.

FeatureList bsfm::PoseEstimator2D3D::points_2d_
private

Definition at line 108 of file pose_estimator_2d3d.h.

Point3DList bsfm::PoseEstimator2D3D::points_3d_
private

Definition at line 109 of file pose_estimator_2d3d.h.

Matrix3d bsfm::PoseEstimator2D3D::T_
private

Definition at line 115 of file pose_estimator_2d3d.h.

Matrix4d bsfm::PoseEstimator2D3D::U_
private

Definition at line 116 of file pose_estimator_2d3d.h.


The documentation for this class was generated from the following files: