Blender V5.0
euclidean_resection.cc
Go to the documentation of this file.
1// Copyright (c) 2009 libmv authors.
2//
3// Permission is hereby granted, free of charge, to any person obtaining a copy
4// of this software and associated documentation files (the "Software"), to
5// deal in the Software without restriction, including without limitation the
6// rights to use, copy, modify, merge, publish, distribute, sublicense, and/or
7// sell copies of the Software, and to permit persons to whom the Software is
8// furnished to do so, subject to the following conditions:
9//
10// The above copyright notice and this permission notice shall be included in
11// all copies or substantial portions of the Software.
12//
13// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
14// IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
15// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
16// AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
17// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
18// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS
19// IN THE SOFTWARE.
20
22
23#include <cmath>
24#include <limits>
25
26#include <Eigen/Geometry>
27#include <Eigen/SVD>
28
29#include "libmv/base/vector.h"
32
33namespace libmv {
35
36typedef unsigned int uint;
37
38bool EuclideanResection(const Mat2X& x_camera,
39 const Mat3X& X_world,
40 Mat3* R,
41 Vec3* t,
42 ResectionMethod method) {
43 switch (method) {
45 EuclideanResectionAnsarDaniilidis(x_camera, X_world, R, t);
46 break;
47 case RESECTION_EPNP:
48 return EuclideanResectionEPnP(x_camera, X_world, R, t);
49 break;
50 case RESECTION_PPNP:
51 return EuclideanResectionPPnP(x_camera, X_world, R, t);
52 break;
53 default: LOG(FATAL) << "Unknown resection method.";
54 }
55 return false;
56}
57
58bool EuclideanResection(const Mat& x_image,
59 const Mat3X& X_world,
60 const Mat3& K,
61 Mat3* R,
62 Vec3* t,
63 ResectionMethod method) {
64 CHECK(x_image.rows() == 2 || x_image.rows() == 3)
65 << "Invalid size for x_image: " << x_image.rows() << "x"
66 << x_image.cols();
67
68 Mat2X x_camera;
69 if (x_image.rows() == 2) {
70 EuclideanToNormalizedCamera(x_image, K, &x_camera);
71 } else if (x_image.rows() == 3) {
72 HomogeneousToNormalizedCamera(x_image, K, &x_camera);
73 }
74 return EuclideanResection(x_camera, X_world, R, t, method);
75}
76
77void AbsoluteOrientation(const Mat3X& X, const Mat3X& Xp, Mat3* R, Vec3* t) {
78 int num_points = X.cols();
79 Vec3 C = X.rowwise().sum() / num_points; // Centroid of X.
80 Vec3 Cp = Xp.rowwise().sum() / num_points; // Centroid of Xp.
81
82 // Normalize the two point sets.
83 Mat3X Xn(3, num_points), Xpn(3, num_points);
84 for (int i = 0; i < num_points; ++i) {
85 Xn.col(i) = X.col(i) - C;
86 Xpn.col(i) = Xp.col(i) - Cp;
87 }
88
89 // Construct the N matrix (pg. 635).
90 double Sxx = Xn.row(0).dot(Xpn.row(0));
91 double Syy = Xn.row(1).dot(Xpn.row(1));
92 double Szz = Xn.row(2).dot(Xpn.row(2));
93 double Sxy = Xn.row(0).dot(Xpn.row(1));
94 double Syx = Xn.row(1).dot(Xpn.row(0));
95 double Sxz = Xn.row(0).dot(Xpn.row(2));
96 double Szx = Xn.row(2).dot(Xpn.row(0));
97 double Syz = Xn.row(1).dot(Xpn.row(2));
98 double Szy = Xn.row(2).dot(Xpn.row(1));
99
100 Mat4 N;
101 // clang-format off
102 N << Sxx + Syy + Szz, Syz - Szy, Szx - Sxz, Sxy - Syx,
103 Syz - Szy, Sxx - Syy - Szz, Sxy + Syx, Szx + Sxz,
104 Szx - Sxz, Sxy + Syx, -Sxx + Syy - Szz, Syz + Szy,
105 Sxy - Syx, Szx + Sxz, Syz + Szy, -Sxx - Syy + Szz;
106 // clang-format on
107
108 // Find the unit quaternion q that maximizes qNq. It is the eigenvector
109 // corresponding to the lagest eigenvalue.
110 Vec4 q = N.jacobiSvd(Eigen::ComputeFullU).matrixU().col(0);
111
112 // Retrieve the 3x3 rotation matrix.
113 Vec4 qq = q.array() * q.array();
114 double q0q1 = q(0) * q(1);
115 double q0q2 = q(0) * q(2);
116 double q0q3 = q(0) * q(3);
117 double q1q2 = q(1) * q(2);
118 double q1q3 = q(1) * q(3);
119 double q2q3 = q(2) * q(3);
120
121 // clang-format off
122 (*R) << qq(0) + qq(1) - qq(2) - qq(3),
123 2 * (q1q2 - q0q3),
124 2 * (q1q3 + q0q2),
125 2 * (q1q2+ q0q3),
126 qq(0) - qq(1) + qq(2) - qq(3),
127 2 * (q2q3 - q0q1),
128 2 * (q1q3 - q0q2),
129 2 * (q2q3 + q0q1),
130 qq(0) - qq(1) - qq(2) + qq(3);
131 // clang-format on
132
133 // Fix the handedness of the R matrix.
134 if (R->determinant() < 0) {
135 R->row(2) = -R->row(2);
136 }
137 // Compute the final translation.
138 *t = Cp - *R * C;
139}
140
141// Convert i and j indices of the original variables into their quadratic
142// permutation single index. It follows that t_ij = t_ji.
143static int IJToPointIndex(int i, int j, int num_points) {
144 // Always make sure that j is bigger than i. This handles t_ij = t_ji.
145 if (j < i) {
146 std::swap(i, j);
147 }
148 int idx;
149 int num_permutation_rows = num_points * (num_points - 1) / 2;
150
151 // All t_ii's are located at the end of the t vector after all t_ij's.
152 if (j == i) {
153 idx = num_permutation_rows + i;
154 } else {
155 int offset = (num_points - i - 1) * (num_points - i) / 2;
156 idx = (num_permutation_rows - offset + j - i - 1);
157 }
158 return idx;
159};
160
161// Convert i and j indexes of the solution for lambda to their linear indexes.
162static int IJToIndex(int i, int j, int num_lambda) {
163 if (j < i) {
164 std::swap(i, j);
165 }
166 int A = num_lambda * (num_lambda + 1) / 2;
167 int B = num_lambda - i;
168 int C = B * (B + 1) / 2;
169 int idx = A - C + j - i;
170 return idx;
171};
172
173static int Sign(double value) {
174 return (value < 0) ? -1 : 1;
175};
176
177// Organizes a square matrix into a single row constraint on the elements of
178// Lambda to create the constraints in equation (5) in "Linear Pose Estimation
179// from Points or Lines", by Ansar, A. and Daniilidis, PAMI 2003. vol. 25, no.
180// 5.
181static Vec MatrixToConstraint(const Mat& A, int num_k_columns, int num_lambda) {
182 Vec C(num_k_columns);
183 C.setZero();
184 int idx = 0;
185 for (int i = 0; i < num_lambda; ++i) {
186 for (int j = i; j < num_lambda; ++j) {
187 C(idx) = A(i, j);
188 if (i != j) {
189 C(idx) += A(j, i);
190 }
191 ++idx;
192 }
193 }
194 return C;
195}
196
197// Normalizes the columns of vectors.
198static void NormalizeColumnVectors(Mat3X* vectors) {
199 int num_columns = vectors->cols();
200 for (int i = 0; i < num_columns; ++i) {
201 vectors->col(i).normalize();
202 }
203}
204
206 const Mat3X& X_world,
207 Mat3* R,
208 Vec3* t) {
209 CHECK(x_camera.cols() == X_world.cols());
210 CHECK(x_camera.cols() > 3);
211
212 int num_points = x_camera.cols();
213
214 // Copy the normalized camera coords into 3 vectors and normalize them so
215 // that they are unit vectors from the camera center.
216 Mat3X x_camera_unit(3, num_points);
217 x_camera_unit.block(0, 0, 2, num_points) = x_camera;
218 x_camera_unit.row(2).setOnes();
219 NormalizeColumnVectors(&x_camera_unit);
220
221 int num_m_rows = num_points * (num_points - 1) / 2;
222 int num_tt_variables = num_points * (num_points + 1) / 2;
223 int num_m_columns = num_tt_variables + 1;
224 Mat M(num_m_columns, num_m_columns);
225 M.setZero();
226 Matu ij_index(num_tt_variables, 2);
227
228 // Create the constraint equations for the t_ij variables (7) and arrange
229 // them into the M matrix (8). Also store the initial (i, j) indices.
230 int row = 0;
231 for (int i = 0; i < num_points; ++i) {
232 for (int j = i + 1; j < num_points; ++j) {
233 M(row, row) = -2 * x_camera_unit.col(i).dot(x_camera_unit.col(j));
234 M(row, num_m_rows + i) = x_camera_unit.col(i).dot(x_camera_unit.col(i));
235 M(row, num_m_rows + j) = x_camera_unit.col(j).dot(x_camera_unit.col(j));
236 Vec3 Xdiff = X_world.col(i) - X_world.col(j);
237 double center_to_point_distance = Xdiff.norm();
238 M(row, num_m_columns - 1) =
239 -center_to_point_distance * center_to_point_distance;
240 ij_index(row, 0) = i;
241 ij_index(row, 1) = j;
242 ++row;
243 }
244 ij_index(i + num_m_rows, 0) = i;
245 ij_index(i + num_m_rows, 1) = i;
246 }
247
248 int num_lambda = num_points + 1; // Dimension of the null space of M.
249 Mat V = M.jacobiSvd(Eigen::ComputeFullV)
250 .matrixV()
251 .block(0, num_m_rows, num_m_columns, num_lambda);
252
253 // TODO(vess): The number of constraint equations in K (num_k_rows) must be
254 // (num_points + 1) * (num_points + 2)/2. This creates a performance issue
255 // for more than 4 points. It is fine for 4 points at the moment with 18
256 // instead of 15 equations.
257 int num_k_rows =
258 num_m_rows +
259 num_points * (num_points * (num_points - 1) / 2 - num_points + 1);
260 int num_k_columns = num_lambda * (num_lambda + 1) / 2;
261 Mat K(num_k_rows, num_k_columns);
262 K.setZero();
263
264 // Construct the first part of the K matrix corresponding to (t_ii, t_jk) for
265 // i != j.
266 int counter_k_row = 0;
267 for (int idx1 = num_m_rows; idx1 < num_tt_variables; ++idx1) {
268 for (int idx2 = 0; idx2 < num_m_rows; ++idx2) {
269 unsigned int i = ij_index(idx1, 0);
270 unsigned int j = ij_index(idx2, 0);
271 unsigned int k = ij_index(idx2, 1);
272
273 if (i != j && i != k) {
274 int idx3 = IJToPointIndex(i, j, num_points);
275 int idx4 = IJToPointIndex(i, k, num_points);
276
277 K.row(counter_k_row) =
278 MatrixToConstraint(V.row(idx1).transpose() * V.row(idx2) -
279 V.row(idx3).transpose() * V.row(idx4),
280 num_k_columns,
281 num_lambda);
282 ++counter_k_row;
283 }
284 }
285 }
286
287 // Construct the second part of the K matrix corresponding to (t_ii,t_jk) for
288 // j==k.
289 for (int idx1 = num_m_rows; idx1 < num_tt_variables; ++idx1) {
290 for (int idx2 = idx1 + 1; idx2 < num_tt_variables; ++idx2) {
291 unsigned int i = ij_index(idx1, 0);
292 unsigned int j = ij_index(idx2, 0);
293 unsigned int k = ij_index(idx2, 1);
294
295 int idx3 = IJToPointIndex(i, j, num_points);
296 int idx4 = IJToPointIndex(i, k, num_points);
297
298 K.row(counter_k_row) =
299 MatrixToConstraint(V.row(idx1).transpose() * V.row(idx2) -
300 V.row(idx3).transpose() * V.row(idx4),
301 num_k_columns,
302 num_lambda);
303 ++counter_k_row;
304 }
305 }
306 Vec L_sq = K.jacobiSvd(Eigen::ComputeFullV).matrixV().col(num_k_columns - 1);
307
308 // Pivot on the largest element for numerical stability. Afterwards recover
309 // the sign of the lambda solution.
310 double max_L_sq_value = fabs(L_sq(IJToIndex(0, 0, num_lambda)));
311 int max_L_sq_index = 1;
312 for (int i = 1; i < num_lambda; ++i) {
313 double abs_sq_value = fabs(L_sq(IJToIndex(i, i, num_lambda)));
314 if (max_L_sq_value < abs_sq_value) {
315 max_L_sq_value = abs_sq_value;
316 max_L_sq_index = i;
317 }
318 }
319 // Ensure positiveness of the largest value corresponding to lambda_ii.
320 L_sq =
321 L_sq * Sign(L_sq(IJToIndex(max_L_sq_index, max_L_sq_index, num_lambda)));
322
323 Vec L(num_lambda);
324 L(max_L_sq_index) =
325 sqrt(L_sq(IJToIndex(max_L_sq_index, max_L_sq_index, num_lambda)));
326
327 for (int i = 0; i < num_lambda; ++i) {
328 if (i != max_L_sq_index) {
329 L(i) = L_sq(IJToIndex(max_L_sq_index, i, num_lambda)) / L(max_L_sq_index);
330 }
331 }
332
333 // Correct the scale using the fact that the last constraint is equal to 1.
334 L = L / (V.row(num_m_columns - 1).dot(L));
335 Vec X = V * L;
336
337 // Recover the distances from the camera center to the 3D points Q.
338 Vec d(num_points);
339 d.setZero();
340 for (int c_point = num_m_rows; c_point < num_tt_variables; ++c_point) {
341 d(c_point - num_m_rows) = sqrt(X(c_point));
342 }
343
344 // Create the 3D points in the camera system.
345 Mat X_cam(3, num_points);
346 for (int c_point = 0; c_point < num_points; ++c_point) {
347 X_cam.col(c_point) = d(c_point) * x_camera_unit.col(c_point);
348 }
349 // Recover the camera translation and rotation.
350 AbsoluteOrientation(X_world, X_cam, R, t);
351}
352
353// Selects 4 virtual control points using mean and PCA.
354static void SelectControlPoints(const Mat3X& X_world,
355 Mat* X_centered,
356 Mat34* X_control_points) {
357 size_t num_points = X_world.cols();
358
359 // The first virtual control point, C0, is the centroid.
360 Vec mean, variance;
361 MeanAndVarianceAlongRows(X_world, &mean, &variance);
362 X_control_points->col(0) = mean;
363
364 // Computes PCA
365 X_centered->resize(3, num_points);
366 for (size_t c = 0; c < num_points; c++) {
367 X_centered->col(c) = X_world.col(c) - mean;
368 }
369 Mat3 X_centered_sq = (*X_centered) * X_centered->transpose();
370 Eigen::JacobiSVD<Mat3> X_centered_sq_svd(X_centered_sq, Eigen::ComputeFullU);
371 Vec3 w = X_centered_sq_svd.singularValues();
372 Mat3 u = X_centered_sq_svd.matrixU();
373 for (size_t c = 0; c < 3; c++) {
374 double k = sqrt(w(c) / num_points);
375 X_control_points->col(c + 1) = mean + k * u.col(c);
376 }
377}
378
379// Computes the barycentric coordinates for all real points
380static bool ComputeBarycentricCoordinates(const Mat3X& X_world_centered,
381 const Mat34& X_control_points,
382 Mat4X* alphas) {
383 size_t num_points = X_world_centered.cols();
384 Mat3 C2;
385 for (size_t c = 1; c < 4; c++) {
386 C2.col(c - 1) = X_control_points.col(c) - X_control_points.col(0);
387 }
388
389 // The selected basis might have correlated vectors causing inverse to produce
390 // undefined result which could lead to issues later on in SVD decomposition.
391 //
392 // Evidently, when built with MSVC it could lead to crash in SVD when the
393 // input contains NaN values:
394 // https://projects.blender.org/blender/blender/issues/145185
395 //
396 // TODO(sergey): Look into choosing a different basis in these cases.
397 bool is_invertible;
398 Mat3 C2inv;
399 C2.computeInverseWithCheck(C2inv, is_invertible);
400 if (!is_invertible) {
401 return false;
402 }
403
404 Mat3X a = C2inv * X_world_centered;
405
406 alphas->resize(4, num_points);
407 alphas->setZero();
408 alphas->block(1, 0, 3, num_points) = a;
409 for (size_t c = 0; c < num_points; c++) {
410 (*alphas)(0, c) = 1.0 - alphas->col(c).sum();
411 }
412
413 return true;
414}
415
416// Estimates the coordinates of all real points in the camera coordinate frame
418 const Mat4X& alphas,
419 const Vec4& betas,
420 const Eigen::Matrix<double, 12, 12>& U,
421 Mat3X* X_camera) {
422 size_t num_points = alphas.cols();
423
424 // Estimates the control points in the camera reference frame.
425 Mat34 C2b;
426 C2b.setZero();
427 for (size_t cu = 0; cu < 4; cu++) {
428 for (size_t c = 0; c < 4; c++) {
429 C2b.col(c) += betas(cu) * U.block(11 - cu, c * 3, 1, 3).transpose();
430 }
431 }
432
433 // Estimates the 3D points in the camera reference frame
434 X_camera->resize(3, num_points);
435 for (size_t c = 0; c < num_points; c++) {
436 X_camera->col(c) = C2b * alphas.col(c);
437 }
438
439 // Check the sign of the z coordinate of the points (should be positive)
440 uint num_z_neg = 0;
441 for (size_t i = 0; i < X_camera->cols(); ++i) {
442 if ((*X_camera)(2, i) < 0) {
443 num_z_neg++;
444 }
445 }
446
447 // If more than 50% of z are negative, we change the signs
448 if (num_z_neg > 0.5 * X_camera->cols()) {
449 C2b = -C2b;
450 *X_camera = -(*X_camera);
451 }
452}
453
454bool EuclideanResectionEPnP(const Mat2X& x_camera,
455 const Mat3X& X_world,
456 Mat3* R,
457 Vec3* t) {
458 CHECK(x_camera.cols() == X_world.cols());
459 CHECK(x_camera.cols() > 3);
460 size_t num_points = X_world.cols();
461
462 // Select the control points.
463 Mat34 X_control_points;
464 Mat X_centered;
465 SelectControlPoints(X_world, &X_centered, &X_control_points);
466
467 // Compute the barycentric coordinates.
468 Mat4X alphas(4, num_points);
469 if (!ComputeBarycentricCoordinates(X_centered, X_control_points, &alphas)) {
470 return false;
471 }
472
473 // Estimates the M matrix with the barycentric coordinates
474 Mat M(2 * num_points, 12);
475 Eigen::Matrix<double, 2, 12> sub_M;
476 for (size_t c = 0; c < num_points; c++) {
477 double a0 = alphas(0, c);
478 double a1 = alphas(1, c);
479 double a2 = alphas(2, c);
480 double a3 = alphas(3, c);
481 double ui = x_camera(0, c);
482 double vi = x_camera(1, c);
483 // clang-format off
484 M.block(2*c, 0, 2, 12) << a0, 0,
485 a0*(-ui), a1, 0,
486 a1*(-ui), a2, 0,
487 a2*(-ui), a3, 0,
488 a3*(-ui), 0,
489 a0, a0*(-vi), 0,
490 a1, a1*(-vi), 0,
491 a2, a2*(-vi), 0,
492 a3, a3*(-vi);
493 // clang-format on
494 }
495
496 // TODO(julien): Avoid the transpose by rewriting the u2.block() calls.
497 Eigen::JacobiSVD<Mat> MtMsvd(M.transpose() * M, Eigen::ComputeFullU);
498 Eigen::Matrix<double, 12, 12> u2 = MtMsvd.matrixU().transpose();
499
500 // Estimate the L matrix.
501 Eigen::Matrix<double, 6, 3> dv1;
502 Eigen::Matrix<double, 6, 3> dv2;
503 Eigen::Matrix<double, 6, 3> dv3;
504 Eigen::Matrix<double, 6, 3> dv4;
505
506 dv1.row(0) = u2.block(11, 0, 1, 3) - u2.block(11, 3, 1, 3);
507 dv1.row(1) = u2.block(11, 0, 1, 3) - u2.block(11, 6, 1, 3);
508 dv1.row(2) = u2.block(11, 0, 1, 3) - u2.block(11, 9, 1, 3);
509 dv1.row(3) = u2.block(11, 3, 1, 3) - u2.block(11, 6, 1, 3);
510 dv1.row(4) = u2.block(11, 3, 1, 3) - u2.block(11, 9, 1, 3);
511 dv1.row(5) = u2.block(11, 6, 1, 3) - u2.block(11, 9, 1, 3);
512 dv2.row(0) = u2.block(10, 0, 1, 3) - u2.block(10, 3, 1, 3);
513 dv2.row(1) = u2.block(10, 0, 1, 3) - u2.block(10, 6, 1, 3);
514 dv2.row(2) = u2.block(10, 0, 1, 3) - u2.block(10, 9, 1, 3);
515 dv2.row(3) = u2.block(10, 3, 1, 3) - u2.block(10, 6, 1, 3);
516 dv2.row(4) = u2.block(10, 3, 1, 3) - u2.block(10, 9, 1, 3);
517 dv2.row(5) = u2.block(10, 6, 1, 3) - u2.block(10, 9, 1, 3);
518 dv3.row(0) = u2.block(9, 0, 1, 3) - u2.block(9, 3, 1, 3);
519 dv3.row(1) = u2.block(9, 0, 1, 3) - u2.block(9, 6, 1, 3);
520 dv3.row(2) = u2.block(9, 0, 1, 3) - u2.block(9, 9, 1, 3);
521 dv3.row(3) = u2.block(9, 3, 1, 3) - u2.block(9, 6, 1, 3);
522 dv3.row(4) = u2.block(9, 3, 1, 3) - u2.block(9, 9, 1, 3);
523 dv3.row(5) = u2.block(9, 6, 1, 3) - u2.block(9, 9, 1, 3);
524 dv4.row(0) = u2.block(8, 0, 1, 3) - u2.block(8, 3, 1, 3);
525 dv4.row(1) = u2.block(8, 0, 1, 3) - u2.block(8, 6, 1, 3);
526 dv4.row(2) = u2.block(8, 0, 1, 3) - u2.block(8, 9, 1, 3);
527 dv4.row(3) = u2.block(8, 3, 1, 3) - u2.block(8, 6, 1, 3);
528 dv4.row(4) = u2.block(8, 3, 1, 3) - u2.block(8, 9, 1, 3);
529 dv4.row(5) = u2.block(8, 6, 1, 3) - u2.block(8, 9, 1, 3);
530
531 Eigen::Matrix<double, 6, 10> L;
532 for (size_t r = 0; r < 6; r++) {
533 // clang-format off
534 L.row(r) << dv1.row(r).dot(dv1.row(r)),
535 2.0 * dv1.row(r).dot(dv2.row(r)),
536 dv2.row(r).dot(dv2.row(r)),
537 2.0 * dv1.row(r).dot(dv3.row(r)),
538 2.0 * dv2.row(r).dot(dv3.row(r)),
539 dv3.row(r).dot(dv3.row(r)),
540 2.0 * dv1.row(r).dot(dv4.row(r)),
541 2.0 * dv2.row(r).dot(dv4.row(r)),
542 2.0 * dv3.row(r).dot(dv4.row(r)),
543 dv4.row(r).dot(dv4.row(r));
544 // clang-format on
545 }
546 Vec6 rho;
547 // clang-format off
548 rho << (X_control_points.col(0) - X_control_points.col(1)).squaredNorm(),
549 (X_control_points.col(0) - X_control_points.col(2)).squaredNorm(),
550 (X_control_points.col(0) - X_control_points.col(3)).squaredNorm(),
551 (X_control_points.col(1) - X_control_points.col(2)).squaredNorm(),
552 (X_control_points.col(1) - X_control_points.col(3)).squaredNorm(),
553 (X_control_points.col(2) - X_control_points.col(3)).squaredNorm();
554 // clang-format on
555
556 // There are three possible solutions based on the three approximations of L
557 // (betas). Below, each one is solved for then the best one is chosen.
558 Mat3X X_camera;
559 Mat3 K;
560 K.setIdentity();
561 vector<Mat3> Rs(3);
562 vector<Vec3> ts(3);
563 Vec rmse(3);
564
565 // At one point this threshold was 1e-3, and caused no end of problems in
566 // Blender by causing frames to not resect when they would have worked fine.
567 // When the resect failed, the projective followup is run leading to worse
568 // results, and often the dreaded "flipping" where objects get flipped
569 // between frames. Instead, disable the check for now, always succeeding. The
570 // ultimate check is always reprojection error anyway.
571 //
572 // TODO(keir): Decide if setting this to infinity, effectively disabling the
573 // check, is the right approach. So far this seems the case.
574 double kSuccessThreshold = std::numeric_limits<double>::max();
575
576 // Find the first possible solution for R, t corresponding to:
577 // Betas = [b00 b01 b11 b02 b12 b22 b03 b13 b23 b33]
578 // Betas_approx_1 = [b00 b01 b02 b03]
579 Vec4 betas = Vec4::Zero();
580 Eigen::Matrix<double, 6, 4> l_6x4;
581 for (size_t r = 0; r < 6; r++) {
582 l_6x4.row(r) << L(r, 0), L(r, 1), L(r, 3), L(r, 6);
583 }
584 Eigen::JacobiSVD<Mat> svd_of_l4(l_6x4,
585 Eigen::ComputeFullU | Eigen::ComputeFullV);
586 Vec4 b4 = svd_of_l4.solve(rho);
587 if ((l_6x4 * b4).isApprox(rho, kSuccessThreshold)) {
588 if (b4(0) < 0) {
589 b4 = -b4;
590 }
591 b4(0) = std::sqrt(b4(0));
592 betas << b4(0), b4(1) / b4(0), b4(2) / b4(0), b4(3) / b4(0);
593 ComputePointsCoordinatesInCameraFrame(alphas, betas, u2, &X_camera);
594 AbsoluteOrientation(X_world, X_camera, &Rs[0], &ts[0]);
595 rmse(0) = RootMeanSquareError(x_camera, X_world, K, Rs[0], ts[0]);
596 } else {
597 LOG(ERROR) << "First approximation of beta not good enough.";
598 ts[0].setZero();
599 rmse(0) = std::numeric_limits<double>::max();
600 }
601
602 // Find the second possible solution for R, t corresponding to:
603 // Betas = [b00 b01 b11 b02 b12 b22 b03 b13 b23 b33]
604 // Betas_approx_2 = [b00 b01 b11]
605 betas.setZero();
606 Eigen::Matrix<double, 6, 3> l_6x3;
607 l_6x3 = L.block(0, 0, 6, 3);
608 Eigen::JacobiSVD<Mat> svdOfL3(l_6x3,
609 Eigen::ComputeFullU | Eigen::ComputeFullV);
610 Vec3 b3 = svdOfL3.solve(rho);
611 VLOG(2) << " rho = " << rho;
612 VLOG(2) << " l_6x3 * b3 = " << l_6x3 * b3;
613 if ((l_6x3 * b3).isApprox(rho, kSuccessThreshold)) {
614 if (b3(0) < 0) {
615 betas(0) = std::sqrt(-b3(0));
616 betas(1) = (b3(2) < 0) ? std::sqrt(-b3(2)) : 0;
617 } else {
618 betas(0) = std::sqrt(b3(0));
619 betas(1) = (b3(2) > 0) ? std::sqrt(b3(2)) : 0;
620 }
621 if (b3(1) < 0) {
622 betas(0) = -betas(0);
623 }
624 betas(2) = 0;
625 betas(3) = 0;
626 ComputePointsCoordinatesInCameraFrame(alphas, betas, u2, &X_camera);
627 AbsoluteOrientation(X_world, X_camera, &Rs[1], &ts[1]);
628 rmse(1) = RootMeanSquareError(x_camera, X_world, K, Rs[1], ts[1]);
629 } else {
630 LOG(ERROR) << "Second approximation of beta not good enough.";
631 ts[1].setZero();
632 rmse(1) = std::numeric_limits<double>::max();
633 }
634
635 // Find the third possible solution for R, t corresponding to:
636 // Betas = [b00 b01 b11 b02 b12 b22 b03 b13 b23 b33]
637 // Betas_approx_3 = [b00 b01 b11 b02 b12]
638 betas.setZero();
639 Eigen::Matrix<double, 6, 5> l_6x5;
640 l_6x5 = L.block(0, 0, 6, 5);
641 Eigen::JacobiSVD<Mat> svdOfL5(l_6x5,
642 Eigen::ComputeFullU | Eigen::ComputeFullV);
643 Vec5 b5 = svdOfL5.solve(rho);
644 if ((l_6x5 * b5).isApprox(rho, kSuccessThreshold)) {
645 if (b5(0) < 0) {
646 betas(0) = std::sqrt(-b5(0));
647 if (b5(2) < 0) {
648 betas(1) = std::sqrt(-b5(2));
649 } else {
650 b5(2) = 0;
651 }
652 } else {
653 betas(0) = std::sqrt(b5(0));
654 if (b5(2) > 0) {
655 betas(1) = std::sqrt(b5(2));
656 } else {
657 b5(2) = 0;
658 }
659 }
660 if (b5(1) < 0) {
661 betas(0) = -betas(0);
662 }
663 betas(2) = b5(3) / betas(0);
664 betas(3) = 0;
665 ComputePointsCoordinatesInCameraFrame(alphas, betas, u2, &X_camera);
666 AbsoluteOrientation(X_world, X_camera, &Rs[2], &ts[2]);
667 rmse(2) = RootMeanSquareError(x_camera, X_world, K, Rs[2], ts[2]);
668 } else {
669 LOG(ERROR) << "Third approximation of beta not good enough.";
670 ts[2].setZero();
671 rmse(2) = std::numeric_limits<double>::max();
672 }
673
674 // Finally, with all three solutions, select the (R, t) with the best RMSE.
675 VLOG(2) << "RMSE for solution 0: " << rmse(0);
676 VLOG(2) << "RMSE for solution 1: " << rmse(1);
677 VLOG(2) << "RMSE for solution 2: " << rmse(2);
678 size_t n = 0;
679 if (rmse(1) < rmse(0)) {
680 n = 1;
681 }
682 if (rmse(2) < rmse(n)) {
683 n = 2;
684 }
685 if (rmse(n) == std::numeric_limits<double>::max()) {
686 LOG(ERROR) << "All three possibilities failed. Reporting failure.";
687 return false;
688 }
689
690 VLOG(1) << "RMSE for best solution #" << n << ": " << rmse(n);
691 *R = Rs[n];
692 *t = ts[n];
693
694 // TODO(julien): Improve the solutions with non-linear refinement.
695 return true;
696}
697
698/*
699
700 Straight from the paper:
701 http://www.diegm.uniud.it/fusiello/papers/3dimpvt12-b.pdf
702
703 function [R T] = ppnp(P,S,tol)
704 % input
705 % P : matrix (nx3) image coordinates in camera reference [u v 1]
706 % S : matrix (nx3) coordinates in world reference [X Y Z]
707 % tol: exit threshold
708 %
709 % output
710 % R : matrix (3x3) rotation (world-to-camera)
711 % T : vector (3x1) translation (world-to-camera)
712 %
713 n = size(P,1);
714 Z = zeros(n);
715 e = ones(n,1);
716 A = eye(n)-((e*e’)./n);
717 II = e./n;
718 err = +Inf;
719 E_old = 1000*ones(n,3);
720 while err>tol
721 [U,˜,V] = svd(P’*Z*A*S);
722 VT = V’;
723 R=U*[1 0 0; 0 1 0; 0 0 det(U*VT)]*VT;
724 PR = P*R;
725 c = (S-Z*PR)’*II;
726 Y = S-e*c’;
727 Zmindiag = diag(PR*Y’)./(sum(P.*P,2));
728 Zmindiag(Zmindiag<0)=0;
729 Z = diag(Zmindiag);
730 E = Y-Z*PR;
731 err = norm(E-E_old,’fro’);
732 E_old = E;
733 end
734 T = -R*c;
735 end
736
737 */
738// TODO(keir): Re-do all the variable names and add comments matching the paper.
739// This implementation has too much of the terseness of the original. On the
740// other hand, it did work on the first try.
741bool EuclideanResectionPPnP(const Mat2X& x_camera,
742 const Mat3X& X_world,
743 Mat3* R,
744 Vec3* t) {
745 int n = x_camera.cols();
746 Mat Z = Mat::Zero(n, n);
747 Vec e = Vec::Ones(n);
748 Mat A = Mat::Identity(n, n) - (e * e.transpose() / n);
749 Vec II = e / n;
750
751 Mat P(n, 3);
752 P.col(0) = x_camera.row(0);
753 P.col(1) = x_camera.row(1);
754 P.col(2).setConstant(1.0);
755
756 Mat S = X_world.transpose();
757
758 double error = std::numeric_limits<double>::infinity();
759 Mat E_old = 1000 * Mat::Ones(n, 3);
760
761 Vec3 c;
762 Mat E(n, 3);
763
764 int iteration = 0;
765 double tolerance = 1e-5;
766 // TODO(keir): The limit of 100 can probably be reduced, but this will require
767 // some investigation.
768 while (error > tolerance && iteration < 100) {
769 Mat3 tmp = P.transpose() * Z * A * S;
770 Eigen::JacobiSVD<Mat3> svd(tmp, Eigen::ComputeFullU | Eigen::ComputeFullV);
771 Mat3 U = svd.matrixU();
772 Mat3 VT = svd.matrixV().transpose();
773 Vec3 s;
774 s << 1, 1, (U * VT).determinant();
775 *R = U * s.asDiagonal() * VT;
776 Mat PR = P * *R; // n x 3
777 c = (S - Z * PR).transpose() * II;
778 Mat Y = S - e * c.transpose(); // n x 3
779 Vec Zmindiag = (PR * Y.transpose())
780 .diagonal()
781 .cwiseQuotient(P.rowwise().squaredNorm());
782 for (int i = 0; i < n; ++i) {
783 Zmindiag[i] = std::max(Zmindiag[i], 0.0);
784 }
785 Z = Zmindiag.asDiagonal();
786 E = Y - Z * PR;
787 error = (E - E_old).norm();
788 LG << "PPnP error(" << (iteration++) << "): " << error;
789 E_old = E;
790 }
791 *t = -*R * c;
792
793 // TODO(keir): Figure out what the failure cases are. Is it too many
794 // iterations? Spend some time going through the math figuring out if there
795 // is some way to detect that the algorithm is going crazy, and return false.
796 return true;
797}
798
799} // namespace euclidean_resection
800} // namespace libmv
#define K(key)
#define X
#define Z
#define Y
#define C
Definition RandGen.cpp:29
#define U
ATTR_WARN_UNUSED_RESULT const BMVert const BMEdge * e
#define A
btScalar determinant() const
Return the determinant of the matrix.
SIMD_FORCE_INLINE const btScalar & w() const
Return the w value.
Definition btQuadWord.h:119
SIMD_FORCE_INLINE btScalar norm() const
Return the norm (length) of the vector.
Definition btVector3.h:263
dot(value.rgb, luminance_coefficients)") DEFINE_VALUE("REDUCE(lhs
#define sqrt
MatBase< R, C > transpose(MatBase< C, R >) RET
#define LG
#define CHECK(expression)
Definition log.h:116
#define LOG(level)
Definition log.h:97
ccl_device_inline float2 fabs(const float2 a)
#define M
#define N
#define B
#define R
#define L
static void error(const char *str)
bool EuclideanResection(const Mat2X &x_camera, const Mat3X &X_world, Mat3 *R, Vec3 *t, ResectionMethod method)
void AbsoluteOrientation(const Mat3X &X, const Mat3X &Xp, Mat3 *R, Vec3 *t)
static void SelectControlPoints(const Mat3X &X_world, Mat *X_centered, Mat34 *X_control_points)
static int IJToIndex(int i, int j, int num_lambda)
static int IJToPointIndex(int i, int j, int num_points)
bool EuclideanResectionEPnP(const Mat2X &x_camera, const Mat3X &X_world, Mat3 *R, Vec3 *t)
static void ComputePointsCoordinatesInCameraFrame(const Mat4X &alphas, const Vec4 &betas, const Eigen::Matrix< double, 12, 12 > &U, Mat3X *X_camera)
static bool ComputeBarycentricCoordinates(const Mat3X &X_world_centered, const Mat34 &X_control_points, Mat4X *alphas)
bool EuclideanResectionPPnP(const Mat2X &x_camera, const Mat3X &X_world, Mat3 *R, Vec3 *t)
static Vec MatrixToConstraint(const Mat &A, int num_k_columns, int num_lambda)
static int Sign(double value)
void EuclideanResectionAnsarDaniilidis(const Mat2X &x_camera, const Mat3X &X_world, Mat3 *R, Vec3 *t)
static void NormalizeColumnVectors(Mat3X *vectors)
Eigen::VectorXd Vec
Definition numeric.h:64
Eigen::Vector4d Vec4
Definition numeric.h:110
Eigen::Matrix< double, 6, 1 > Vec6
Definition numeric.h:112
Eigen::Matrix< double, 4, 4 > Mat4
Definition numeric.h:80
Eigen::Matrix< double, 3, 3 > Mat3
Definition numeric.h:75
void HomogeneousToNormalizedCamera(const Mat3X &x, const Mat3 &K, Mat2X *n)
Eigen::MatrixXd Mat
Definition numeric.h:63
void MeanAndVarianceAlongRows(const Mat &A, Vec *mean_pointer, Vec *variance_pointer)
Definition numeric.cc:90
double RootMeanSquareError(const Mat2X &x_image, const Mat4X &X_world, const Mat34 &P)
Estimates the root mean square error (2D).
Eigen::Matrix< double, 3, 4 > Mat34
Definition numeric.h:76
std::vector< ElementType, Eigen::aligned_allocator< ElementType > > vector
Eigen::Matrix< double, 5, 1 > Vec5
Definition numeric.h:111
Eigen::Vector3d Vec3
Definition numeric.h:109
Eigen::Matrix< double, 4, Eigen::Dynamic > Mat4X
Definition numeric.h:96
Eigen::Matrix< double, 3, Eigen::Dynamic > Mat3X
Definition numeric.h:95
Eigen::Matrix< unsigned int, Eigen::Dynamic, Eigen::Dynamic > Matu
Definition numeric.h:69
Eigen::Matrix< double, 2, Eigen::Dynamic > Mat2X
Definition numeric.h:94
void EuclideanToNormalizedCamera(const Mat2X &x, const Mat3 &K, Mat2X *n)
i
Definition text_draw.cc:230
CCL_NAMESPACE_BEGIN struct Window V