blob: 18b7e8c0f8ccfa5684da4a65a3556460297ffb04 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2015 Google Inc. All rights reserved.
3// http://ceres-solver.org/
4//
5// Redistribution and use in source and binary forms, with or without
6// modification, are permitted provided that the following conditions are met:
7//
8// * Redistributions of source code must retain the above copyright notice,
9// this list of conditions and the following disclaimer.
10// * Redistributions in binary form must reproduce the above copyright notice,
11// this list of conditions and the following disclaimer in the documentation
12// and/or other materials provided with the distribution.
13// * Neither the name of Google Inc. nor the names of its contributors may be
14// used to endorse or promote products derived from this software without
15// specific prior written permission.
16//
17// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20// ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21// LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22// CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23// SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24// INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25// CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26// ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27// POSSIBILITY OF SUCH DAMAGE.
28//
29// Author: sameeragarwal@google.com (Sameer Agarwal)
30
31#include <cmath>
32#include <limits>
33#include <memory>
34
35#include "Eigen/Geometry"
36#include "ceres/autodiff_local_parameterization.h"
37#include "ceres/householder_vector.h"
38#include "ceres/internal/autodiff.h"
39#include "ceres/internal/eigen.h"
40#include "ceres/local_parameterization.h"
41#include "ceres/random.h"
42#include "ceres/rotation.h"
43#include "gtest/gtest.h"
44
45namespace ceres {
46namespace internal {
47
48TEST(IdentityParameterization, EverythingTest) {
49 IdentityParameterization parameterization(3);
50 EXPECT_EQ(parameterization.GlobalSize(), 3);
51 EXPECT_EQ(parameterization.LocalSize(), 3);
52
53 double x[3] = {1.0, 2.0, 3.0};
54 double delta[3] = {0.0, 1.0, 2.0};
55 double x_plus_delta[3] = {0.0, 0.0, 0.0};
56 parameterization.Plus(x, delta, x_plus_delta);
57 EXPECT_EQ(x_plus_delta[0], 1.0);
58 EXPECT_EQ(x_plus_delta[1], 3.0);
59 EXPECT_EQ(x_plus_delta[2], 5.0);
60
61 double jacobian[9];
62 parameterization.ComputeJacobian(x, jacobian);
63 int k = 0;
64 for (int i = 0; i < 3; ++i) {
65 for (int j = 0; j < 3; ++j, ++k) {
66 EXPECT_EQ(jacobian[k], (i == j) ? 1.0 : 0.0);
67 }
68 }
69
70 Matrix global_matrix = Matrix::Ones(10, 3);
71 Matrix local_matrix = Matrix::Zero(10, 3);
72 parameterization.MultiplyByJacobian(x,
73 10,
74 global_matrix.data(),
75 local_matrix.data());
76 EXPECT_EQ((local_matrix - global_matrix).norm(), 0.0);
77}
78
79
80TEST(SubsetParameterization, NegativeParameterIndexDeathTest) {
81 std::vector<int> constant_parameters;
82 constant_parameters.push_back(-1);
83 EXPECT_DEATH_IF_SUPPORTED(
84 SubsetParameterization parameterization(2, constant_parameters),
85 "greater than equal to zero");
86}
87
88TEST(SubsetParameterization, GreaterThanSizeParameterIndexDeathTest) {
89 std::vector<int> constant_parameters;
90 constant_parameters.push_back(2);
91 EXPECT_DEATH_IF_SUPPORTED(
92 SubsetParameterization parameterization(2, constant_parameters),
93 "less than the size");
94}
95
96TEST(SubsetParameterization, DuplicateParametersDeathTest) {
97 std::vector<int> constant_parameters;
98 constant_parameters.push_back(1);
99 constant_parameters.push_back(1);
100 EXPECT_DEATH_IF_SUPPORTED(
101 SubsetParameterization parameterization(2, constant_parameters),
102 "duplicates");
103}
104
105TEST(SubsetParameterization,
106 ProductParameterizationWithZeroLocalSizeSubsetParameterization1) {
107 std::vector<int> constant_parameters;
108 constant_parameters.push_back(0);
109 LocalParameterization* subset_param =
110 new SubsetParameterization(1, constant_parameters);
111 LocalParameterization* identity_param = new IdentityParameterization(2);
112 ProductParameterization product_param(subset_param, identity_param);
113 EXPECT_EQ(product_param.GlobalSize(), 3);
114 EXPECT_EQ(product_param.LocalSize(), 2);
115 double x[] = {1.0, 1.0, 1.0};
116 double delta[] = {2.0, 3.0};
117 double x_plus_delta[] = {0.0, 0.0, 0.0};
118 EXPECT_TRUE(product_param.Plus(x, delta, x_plus_delta));
119 EXPECT_EQ(x_plus_delta[0], x[0]);
120 EXPECT_EQ(x_plus_delta[1], x[1] + delta[0]);
121 EXPECT_EQ(x_plus_delta[2], x[2] + delta[1]);
122
123 Matrix actual_jacobian(3, 2);
124 EXPECT_TRUE(product_param.ComputeJacobian(x, actual_jacobian.data()));
125}
126
127TEST(SubsetParameterization,
128 ProductParameterizationWithZeroLocalSizeSubsetParameterization2) {
129 std::vector<int> constant_parameters;
130 constant_parameters.push_back(0);
131 LocalParameterization* subset_param =
132 new SubsetParameterization(1, constant_parameters);
133 LocalParameterization* identity_param = new IdentityParameterization(2);
134 ProductParameterization product_param(identity_param, subset_param);
135 EXPECT_EQ(product_param.GlobalSize(), 3);
136 EXPECT_EQ(product_param.LocalSize(), 2);
137 double x[] = {1.0, 1.0, 1.0};
138 double delta[] = {2.0, 3.0};
139 double x_plus_delta[] = {0.0, 0.0, 0.0};
140 EXPECT_TRUE(product_param.Plus(x, delta, x_plus_delta));
141 EXPECT_EQ(x_plus_delta[0], x[0] + delta[0]);
142 EXPECT_EQ(x_plus_delta[1], x[1] + delta[1]);
143 EXPECT_EQ(x_plus_delta[2], x[2]);
144
145 Matrix actual_jacobian(3, 2);
146 EXPECT_TRUE(product_param.ComputeJacobian(x, actual_jacobian.data()));
147}
148
149TEST(SubsetParameterization, NormalFunctionTest) {
150 const int kGlobalSize = 4;
151 const int kLocalSize = 3;
152
153 double x[kGlobalSize] = {1.0, 2.0, 3.0, 4.0};
154 for (int i = 0; i < kGlobalSize; ++i) {
155 std::vector<int> constant_parameters;
156 constant_parameters.push_back(i);
157 SubsetParameterization parameterization(kGlobalSize, constant_parameters);
158 double delta[kLocalSize] = {1.0, 2.0, 3.0};
159 double x_plus_delta[kGlobalSize] = {0.0, 0.0, 0.0};
160
161 parameterization.Plus(x, delta, x_plus_delta);
162 int k = 0;
163 for (int j = 0; j < kGlobalSize; ++j) {
164 if (j == i) {
165 EXPECT_EQ(x_plus_delta[j], x[j]);
166 } else {
167 EXPECT_EQ(x_plus_delta[j], x[j] + delta[k++]);
168 }
169 }
170
171 double jacobian[kGlobalSize * kLocalSize];
172 parameterization.ComputeJacobian(x, jacobian);
173 int delta_cursor = 0;
174 int jacobian_cursor = 0;
175 for (int j = 0; j < kGlobalSize; ++j) {
176 if (j != i) {
177 for (int k = 0; k < kLocalSize; ++k, jacobian_cursor++) {
178 EXPECT_EQ(jacobian[jacobian_cursor], delta_cursor == k ? 1.0 : 0.0);
179 }
180 ++delta_cursor;
181 } else {
182 for (int k = 0; k < kLocalSize; ++k, jacobian_cursor++) {
183 EXPECT_EQ(jacobian[jacobian_cursor], 0.0);
184 }
185 }
186 }
187
188 Matrix global_matrix = Matrix::Ones(10, kGlobalSize);
189 for (int row = 0; row < kGlobalSize; ++row) {
190 for (int col = 0; col < kGlobalSize; ++col) {
191 global_matrix(row, col) = col;
192 }
193 }
194
195 Matrix local_matrix = Matrix::Zero(10, kLocalSize);
196 parameterization.MultiplyByJacobian(x,
197 10,
198 global_matrix.data(),
199 local_matrix.data());
200 Matrix expected_local_matrix =
201 global_matrix * MatrixRef(jacobian, kGlobalSize, kLocalSize);
202 EXPECT_EQ((local_matrix - expected_local_matrix).norm(), 0.0);
203 }
204}
205
206// Functor needed to implement automatically differentiated Plus for
207// quaternions.
208struct QuaternionPlus {
209 template<typename T>
210 bool operator()(const T* x, const T* delta, T* x_plus_delta) const {
211 const T squared_norm_delta =
212 delta[0] * delta[0] + delta[1] * delta[1] + delta[2] * delta[2];
213
214 T q_delta[4];
215 if (squared_norm_delta > T(0.0)) {
216 T norm_delta = sqrt(squared_norm_delta);
217 const T sin_delta_by_delta = sin(norm_delta) / norm_delta;
218 q_delta[0] = cos(norm_delta);
219 q_delta[1] = sin_delta_by_delta * delta[0];
220 q_delta[2] = sin_delta_by_delta * delta[1];
221 q_delta[3] = sin_delta_by_delta * delta[2];
222 } else {
223 // We do not just use q_delta = [1,0,0,0] here because that is a
224 // constant and when used for automatic differentiation will
225 // lead to a zero derivative. Instead we take a first order
226 // approximation and evaluate it at zero.
227 q_delta[0] = T(1.0);
228 q_delta[1] = delta[0];
229 q_delta[2] = delta[1];
230 q_delta[3] = delta[2];
231 }
232
233 QuaternionProduct(q_delta, x, x_plus_delta);
234 return true;
235 }
236};
237
238template<typename Parameterization, typename Plus>
239void QuaternionParameterizationTestHelper(
240 const double* x, const double* delta,
241 const double* x_plus_delta_ref) {
242 const int kGlobalSize = 4;
243 const int kLocalSize = 3;
244
245 const double kTolerance = 1e-14;
246
247 double x_plus_delta[kGlobalSize] = {0.0, 0.0, 0.0, 0.0};
248 Parameterization parameterization;
249 parameterization.Plus(x, delta, x_plus_delta);
250 for (int i = 0; i < kGlobalSize; ++i) {
251 EXPECT_NEAR(x_plus_delta[i], x_plus_delta[i], kTolerance);
252 }
253
254 const double x_plus_delta_norm =
255 sqrt(x_plus_delta[0] * x_plus_delta[0] +
256 x_plus_delta[1] * x_plus_delta[1] +
257 x_plus_delta[2] * x_plus_delta[2] +
258 x_plus_delta[3] * x_plus_delta[3]);
259
260 EXPECT_NEAR(x_plus_delta_norm, 1.0, kTolerance);
261
262 double jacobian_ref[12];
263 double zero_delta[kLocalSize] = {0.0, 0.0, 0.0};
264 const double* parameters[2] = {x, zero_delta};
265 double* jacobian_array[2] = { NULL, jacobian_ref };
266
267 // Autodiff jacobian at delta_x = 0.
268 internal::AutoDifferentiate<StaticParameterDims<kGlobalSize, kLocalSize>>(
269 Plus(),
270 parameters,
271 kGlobalSize,
272 x_plus_delta,
273 jacobian_array);
274
275 double jacobian[12];
276 parameterization.ComputeJacobian(x, jacobian);
277 for (int i = 0; i < 12; ++i) {
278 EXPECT_TRUE(IsFinite(jacobian[i]));
279 EXPECT_NEAR(jacobian[i], jacobian_ref[i], kTolerance)
280 << "Jacobian mismatch: i = " << i
281 << "\n Expected \n"
282 << ConstMatrixRef(jacobian_ref, kGlobalSize, kLocalSize)
283 << "\n Actual \n"
284 << ConstMatrixRef(jacobian, kGlobalSize, kLocalSize);
285 }
286
287 Matrix global_matrix = Matrix::Random(10, kGlobalSize);
288 Matrix local_matrix = Matrix::Zero(10, kLocalSize);
289 parameterization.MultiplyByJacobian(x,
290 10,
291 global_matrix.data(),
292 local_matrix.data());
293 Matrix expected_local_matrix =
294 global_matrix * MatrixRef(jacobian, kGlobalSize, kLocalSize);
295 EXPECT_NEAR((local_matrix - expected_local_matrix).norm(),
296 0.0,
297 10.0 * std::numeric_limits<double>::epsilon());
298}
299
300template <int N>
301void Normalize(double* x) {
302 VectorRef(x, N).normalize();
303}
304
305TEST(QuaternionParameterization, ZeroTest) {
306 double x[4] = {0.5, 0.5, 0.5, 0.5};
307 double delta[3] = {0.0, 0.0, 0.0};
308 double q_delta[4] = {1.0, 0.0, 0.0, 0.0};
309 double x_plus_delta[4] = {0.0, 0.0, 0.0, 0.0};
310 QuaternionProduct(q_delta, x, x_plus_delta);
311 QuaternionParameterizationTestHelper<QuaternionParameterization,
312 QuaternionPlus>(x, delta, x_plus_delta);
313}
314
315TEST(QuaternionParameterization, NearZeroTest) {
316 double x[4] = {0.52, 0.25, 0.15, 0.45};
317 Normalize<4>(x);
318
319 double delta[3] = {0.24, 0.15, 0.10};
320 for (int i = 0; i < 3; ++i) {
321 delta[i] = delta[i] * 1e-14;
322 }
323
324 double q_delta[4];
325 q_delta[0] = 1.0;
326 q_delta[1] = delta[0];
327 q_delta[2] = delta[1];
328 q_delta[3] = delta[2];
329
330 double x_plus_delta[4] = {0.0, 0.0, 0.0, 0.0};
331 QuaternionProduct(q_delta, x, x_plus_delta);
332 QuaternionParameterizationTestHelper<QuaternionParameterization,
333 QuaternionPlus>(x, delta, x_plus_delta);
334}
335
336TEST(QuaternionParameterization, AwayFromZeroTest) {
337 double x[4] = {0.52, 0.25, 0.15, 0.45};
338 Normalize<4>(x);
339
340 double delta[3] = {0.24, 0.15, 0.10};
341 const double delta_norm = sqrt(delta[0] * delta[0] +
342 delta[1] * delta[1] +
343 delta[2] * delta[2]);
344 double q_delta[4];
345 q_delta[0] = cos(delta_norm);
346 q_delta[1] = sin(delta_norm) / delta_norm * delta[0];
347 q_delta[2] = sin(delta_norm) / delta_norm * delta[1];
348 q_delta[3] = sin(delta_norm) / delta_norm * delta[2];
349
350 double x_plus_delta[4] = {0.0, 0.0, 0.0, 0.0};
351 QuaternionProduct(q_delta, x, x_plus_delta);
352 QuaternionParameterizationTestHelper<QuaternionParameterization,
353 QuaternionPlus>(x, delta, x_plus_delta);
354}
355
356// Functor needed to implement automatically differentiated Plus for
357// Eigen's quaternion.
358struct EigenQuaternionPlus {
359 template<typename T>
360 bool operator()(const T* x, const T* delta, T* x_plus_delta) const {
361 const T norm_delta =
362 sqrt(delta[0] * delta[0] + delta[1] * delta[1] + delta[2] * delta[2]);
363
364 Eigen::Quaternion<T> q_delta;
365 if (norm_delta > T(0.0)) {
366 const T sin_delta_by_delta = sin(norm_delta) / norm_delta;
367 q_delta.coeffs() << sin_delta_by_delta * delta[0],
368 sin_delta_by_delta * delta[1], sin_delta_by_delta * delta[2],
369 cos(norm_delta);
370 } else {
371 // We do not just use q_delta = [0,0,0,1] here because that is a
372 // constant and when used for automatic differentiation will
373 // lead to a zero derivative. Instead we take a first order
374 // approximation and evaluate it at zero.
375 q_delta.coeffs() << delta[0], delta[1], delta[2], T(1.0);
376 }
377
378 Eigen::Map<Eigen::Quaternion<T>> x_plus_delta_ref(x_plus_delta);
379 Eigen::Map<const Eigen::Quaternion<T>> x_ref(x);
380 x_plus_delta_ref = q_delta * x_ref;
381 return true;
382 }
383};
384
385TEST(EigenQuaternionParameterization, ZeroTest) {
386 Eigen::Quaterniond x(0.5, 0.5, 0.5, 0.5);
387 double delta[3] = {0.0, 0.0, 0.0};
388 Eigen::Quaterniond q_delta(1.0, 0.0, 0.0, 0.0);
389 Eigen::Quaterniond x_plus_delta = q_delta * x;
390 QuaternionParameterizationTestHelper<EigenQuaternionParameterization,
391 EigenQuaternionPlus>(
392 x.coeffs().data(), delta, x_plus_delta.coeffs().data());
393}
394
395TEST(EigenQuaternionParameterization, NearZeroTest) {
396 Eigen::Quaterniond x(0.52, 0.25, 0.15, 0.45);
397 x.normalize();
398
399 double delta[3] = {0.24, 0.15, 0.10};
400 for (int i = 0; i < 3; ++i) {
401 delta[i] = delta[i] * 1e-14;
402 }
403
404 // Note: w is first in the constructor.
405 Eigen::Quaterniond q_delta(1.0, delta[0], delta[1], delta[2]);
406
407 Eigen::Quaterniond x_plus_delta = q_delta * x;
408 QuaternionParameterizationTestHelper<EigenQuaternionParameterization,
409 EigenQuaternionPlus>(
410 x.coeffs().data(), delta, x_plus_delta.coeffs().data());
411}
412
413TEST(EigenQuaternionParameterization, AwayFromZeroTest) {
414 Eigen::Quaterniond x(0.52, 0.25, 0.15, 0.45);
415 x.normalize();
416
417 double delta[3] = {0.24, 0.15, 0.10};
418 const double delta_norm = sqrt(delta[0] * delta[0] +
419 delta[1] * delta[1] +
420 delta[2] * delta[2]);
421
422 // Note: w is first in the constructor.
423 Eigen::Quaterniond q_delta(cos(delta_norm),
424 sin(delta_norm) / delta_norm * delta[0],
425 sin(delta_norm) / delta_norm * delta[1],
426 sin(delta_norm) / delta_norm * delta[2]);
427
428 Eigen::Quaterniond x_plus_delta = q_delta * x;
429 QuaternionParameterizationTestHelper<EigenQuaternionParameterization,
430 EigenQuaternionPlus>(
431 x.coeffs().data(), delta, x_plus_delta.coeffs().data());
432}
433
434// Functor needed to implement automatically differentiated Plus for
435// homogeneous vectors. Note this explicitly defined for vectors of size 4.
436struct HomogeneousVectorParameterizationPlus {
437 template<typename Scalar>
438 bool operator()(const Scalar* p_x, const Scalar* p_delta,
439 Scalar* p_x_plus_delta) const {
440 Eigen::Map<const Eigen::Matrix<Scalar, 4, 1>> x(p_x);
441 Eigen::Map<const Eigen::Matrix<Scalar, 3, 1>> delta(p_delta);
442 Eigen::Map<Eigen::Matrix<Scalar, 4, 1>> x_plus_delta(p_x_plus_delta);
443
444 const Scalar squared_norm_delta =
445 delta[0] * delta[0] + delta[1] * delta[1] + delta[2] * delta[2];
446
447 Eigen::Matrix<Scalar, 4, 1> y;
448 Scalar one_half(0.5);
449 if (squared_norm_delta > Scalar(0.0)) {
450 Scalar norm_delta = sqrt(squared_norm_delta);
451 Scalar norm_delta_div_2 = 0.5 * norm_delta;
452 const Scalar sin_delta_by_delta = sin(norm_delta_div_2) /
453 norm_delta_div_2;
454 y[0] = sin_delta_by_delta * delta[0] * one_half;
455 y[1] = sin_delta_by_delta * delta[1] * one_half;
456 y[2] = sin_delta_by_delta * delta[2] * one_half;
457 y[3] = cos(norm_delta_div_2);
458
459 } else {
460 // We do not just use y = [0,0,0,1] here because that is a
461 // constant and when used for automatic differentiation will
462 // lead to a zero derivative. Instead we take a first order
463 // approximation and evaluate it at zero.
464 y[0] = delta[0] * one_half;
465 y[1] = delta[1] * one_half;
466 y[2] = delta[2] * one_half;
467 y[3] = Scalar(1.0);
468 }
469
470 Eigen::Matrix<Scalar, Eigen::Dynamic, 1> v(4);
471 Scalar beta;
472 internal::ComputeHouseholderVector<Scalar>(x, &v, &beta);
473
474 x_plus_delta = x.norm() * (y - v * (beta * v.dot(y)));
475
476 return true;
477 }
478};
479
480void HomogeneousVectorParameterizationHelper(const double* x,
481 const double* delta) {
482 const double kTolerance = 1e-14;
483
484 HomogeneousVectorParameterization homogeneous_vector_parameterization(4);
485
486 // Ensure the update maintains the norm.
487 double x_plus_delta[4] = {0.0, 0.0, 0.0, 0.0};
488 homogeneous_vector_parameterization.Plus(x, delta, x_plus_delta);
489
490 const double x_plus_delta_norm =
491 sqrt(x_plus_delta[0] * x_plus_delta[0] +
492 x_plus_delta[1] * x_plus_delta[1] +
493 x_plus_delta[2] * x_plus_delta[2] +
494 x_plus_delta[3] * x_plus_delta[3]);
495
496 const double x_norm = sqrt(x[0] * x[0] + x[1] * x[1] +
497 x[2] * x[2] + x[3] * x[3]);
498
499 EXPECT_NEAR(x_plus_delta_norm, x_norm, kTolerance);
500
501 // Autodiff jacobian at delta_x = 0.
502 AutoDiffLocalParameterization<HomogeneousVectorParameterizationPlus, 4, 3>
503 autodiff_jacobian;
504
505 double jacobian_autodiff[12];
506 double jacobian_analytic[12];
507
508 homogeneous_vector_parameterization.ComputeJacobian(x, jacobian_analytic);
509 autodiff_jacobian.ComputeJacobian(x, jacobian_autodiff);
510
511 for (int i = 0; i < 12; ++i) {
512 EXPECT_TRUE(ceres::IsFinite(jacobian_analytic[i]));
513 EXPECT_NEAR(jacobian_analytic[i], jacobian_autodiff[i], kTolerance)
514 << "Jacobian mismatch: i = " << i << ", " << jacobian_analytic[i] << " "
515 << jacobian_autodiff[i];
516 }
517}
518
519TEST(HomogeneousVectorParameterization, ZeroTest) {
520 double x[4] = {0.0, 0.0, 0.0, 1.0};
521 Normalize<4>(x);
522 double delta[3] = {0.0, 0.0, 0.0};
523
524 HomogeneousVectorParameterizationHelper(x, delta);
525}
526
527TEST(HomogeneousVectorParameterization, NearZeroTest1) {
528 double x[4] = {1e-5, 1e-5, 1e-5, 1.0};
529 Normalize<4>(x);
530 double delta[3] = {0.0, 1.0, 0.0};
531
532 HomogeneousVectorParameterizationHelper(x, delta);
533}
534
535TEST(HomogeneousVectorParameterization, NearZeroTest2) {
536 double x[4] = {0.001, 0.0, 0.0, 0.0};
537 double delta[3] = {0.0, 1.0, 0.0};
538
539 HomogeneousVectorParameterizationHelper(x, delta);
540}
541
542TEST(HomogeneousVectorParameterization, AwayFromZeroTest1) {
543 double x[4] = {0.52, 0.25, 0.15, 0.45};
544 Normalize<4>(x);
545 double delta[3] = {0.0, 1.0, -0.5};
546
547 HomogeneousVectorParameterizationHelper(x, delta);
548}
549
550TEST(HomogeneousVectorParameterization, AwayFromZeroTest2) {
551 double x[4] = {0.87, -0.25, -0.34, 0.45};
552 Normalize<4>(x);
553 double delta[3] = {0.0, 0.0, -0.5};
554
555 HomogeneousVectorParameterizationHelper(x, delta);
556}
557
558TEST(HomogeneousVectorParameterization, AwayFromZeroTest3) {
559 double x[4] = {0.0, 0.0, 0.0, 2.0};
560 double delta[3] = {0.0, 0.0, 0};
561
562 HomogeneousVectorParameterizationHelper(x, delta);
563}
564
565TEST(HomogeneousVectorParameterization, AwayFromZeroTest4) {
566 double x[4] = {0.2, -1.0, 0.0, 2.0};
567 double delta[3] = {1.4, 0.0, -0.5};
568
569 HomogeneousVectorParameterizationHelper(x, delta);
570}
571
572TEST(HomogeneousVectorParameterization, AwayFromZeroTest5) {
573 double x[4] = {2.0, 0.0, 0.0, 0.0};
574 double delta[3] = {1.4, 0.0, -0.5};
575
576 HomogeneousVectorParameterizationHelper(x, delta);
577}
578
579TEST(HomogeneousVectorParameterization, DeathTests) {
580 EXPECT_DEATH_IF_SUPPORTED(HomogeneousVectorParameterization x(1), "size");
581}
582
583
584class ProductParameterizationTest : public ::testing::Test {
585 protected :
586 virtual void SetUp() {
587 const int global_size1 = 5;
588 std::vector<int> constant_parameters1;
589 constant_parameters1.push_back(2);
590 param1_.reset(new SubsetParameterization(global_size1,
591 constant_parameters1));
592
593 const int global_size2 = 3;
594 std::vector<int> constant_parameters2;
595 constant_parameters2.push_back(0);
596 constant_parameters2.push_back(1);
597 param2_.reset(new SubsetParameterization(global_size2,
598 constant_parameters2));
599
600 const int global_size3 = 4;
601 std::vector<int> constant_parameters3;
602 constant_parameters3.push_back(1);
603 param3_.reset(new SubsetParameterization(global_size3,
604 constant_parameters3));
605
606 const int global_size4 = 2;
607 std::vector<int> constant_parameters4;
608 constant_parameters4.push_back(1);
609 param4_.reset(new SubsetParameterization(global_size4,
610 constant_parameters4));
611 }
612
613 std::unique_ptr<LocalParameterization> param1_;
614 std::unique_ptr<LocalParameterization> param2_;
615 std::unique_ptr<LocalParameterization> param3_;
616 std::unique_ptr<LocalParameterization> param4_;
617};
618
619TEST_F(ProductParameterizationTest, LocalAndGlobalSize2) {
620 LocalParameterization* param1 = param1_.release();
621 LocalParameterization* param2 = param2_.release();
622
623 ProductParameterization product_param(param1, param2);
624 EXPECT_EQ(product_param.LocalSize(),
625 param1->LocalSize() + param2->LocalSize());
626 EXPECT_EQ(product_param.GlobalSize(),
627 param1->GlobalSize() + param2->GlobalSize());
628}
629
630
631TEST_F(ProductParameterizationTest, LocalAndGlobalSize3) {
632 LocalParameterization* param1 = param1_.release();
633 LocalParameterization* param2 = param2_.release();
634 LocalParameterization* param3 = param3_.release();
635
636 ProductParameterization product_param(param1, param2, param3);
637 EXPECT_EQ(product_param.LocalSize(),
638 param1->LocalSize() + param2->LocalSize() + param3->LocalSize());
639 EXPECT_EQ(product_param.GlobalSize(),
640 param1->GlobalSize() + param2->GlobalSize() + param3->GlobalSize());
641}
642
643TEST_F(ProductParameterizationTest, LocalAndGlobalSize4) {
644 LocalParameterization* param1 = param1_.release();
645 LocalParameterization* param2 = param2_.release();
646 LocalParameterization* param3 = param3_.release();
647 LocalParameterization* param4 = param4_.release();
648
649 ProductParameterization product_param(param1, param2, param3, param4);
650 EXPECT_EQ(product_param.LocalSize(),
651 param1->LocalSize() +
652 param2->LocalSize() +
653 param3->LocalSize() +
654 param4->LocalSize());
655 EXPECT_EQ(product_param.GlobalSize(),
656 param1->GlobalSize() +
657 param2->GlobalSize() +
658 param3->GlobalSize() +
659 param4->GlobalSize());
660}
661
662TEST_F(ProductParameterizationTest, Plus) {
663 LocalParameterization* param1 = param1_.release();
664 LocalParameterization* param2 = param2_.release();
665 LocalParameterization* param3 = param3_.release();
666 LocalParameterization* param4 = param4_.release();
667
668 ProductParameterization product_param(param1, param2, param3, param4);
669 std::vector<double> x(product_param.GlobalSize(), 0.0);
670 std::vector<double> delta(product_param.LocalSize(), 0.0);
671 std::vector<double> x_plus_delta_expected(product_param.GlobalSize(), 0.0);
672 std::vector<double> x_plus_delta(product_param.GlobalSize(), 0.0);
673
674 for (int i = 0; i < product_param.GlobalSize(); ++i) {
675 x[i] = RandNormal();
676 }
677
678 for (int i = 0; i < product_param.LocalSize(); ++i) {
679 delta[i] = RandNormal();
680 }
681
682 EXPECT_TRUE(product_param.Plus(&x[0], &delta[0], &x_plus_delta_expected[0]));
683 int x_cursor = 0;
684 int delta_cursor = 0;
685
686 EXPECT_TRUE(param1->Plus(&x[x_cursor],
687 &delta[delta_cursor],
688 &x_plus_delta[x_cursor]));
689 x_cursor += param1->GlobalSize();
690 delta_cursor += param1->LocalSize();
691
692 EXPECT_TRUE(param2->Plus(&x[x_cursor],
693 &delta[delta_cursor],
694 &x_plus_delta[x_cursor]));
695 x_cursor += param2->GlobalSize();
696 delta_cursor += param2->LocalSize();
697
698 EXPECT_TRUE(param3->Plus(&x[x_cursor],
699 &delta[delta_cursor],
700 &x_plus_delta[x_cursor]));
701 x_cursor += param3->GlobalSize();
702 delta_cursor += param3->LocalSize();
703
704 EXPECT_TRUE(param4->Plus(&x[x_cursor],
705 &delta[delta_cursor],
706 &x_plus_delta[x_cursor]));
707 x_cursor += param4->GlobalSize();
708 delta_cursor += param4->LocalSize();
709
710 for (int i = 0; i < x.size(); ++i) {
711 EXPECT_EQ(x_plus_delta[i], x_plus_delta_expected[i]);
712 }
713}
714
715TEST_F(ProductParameterizationTest, ComputeJacobian) {
716 LocalParameterization* param1 = param1_.release();
717 LocalParameterization* param2 = param2_.release();
718 LocalParameterization* param3 = param3_.release();
719 LocalParameterization* param4 = param4_.release();
720
721 ProductParameterization product_param(param1, param2, param3, param4);
722 std::vector<double> x(product_param.GlobalSize(), 0.0);
723
724 for (int i = 0; i < product_param.GlobalSize(); ++i) {
725 x[i] = RandNormal();
726 }
727
728 Matrix jacobian = Matrix::Random(product_param.GlobalSize(),
729 product_param.LocalSize());
730 EXPECT_TRUE(product_param.ComputeJacobian(&x[0], jacobian.data()));
731 int x_cursor = 0;
732 int delta_cursor = 0;
733
734 Matrix jacobian1(param1->GlobalSize(), param1->LocalSize());
735 EXPECT_TRUE(param1->ComputeJacobian(&x[x_cursor], jacobian1.data()));
736 jacobian.block(x_cursor, delta_cursor,
737 param1->GlobalSize(),
738 param1->LocalSize())
739 -= jacobian1;
740 x_cursor += param1->GlobalSize();
741 delta_cursor += param1->LocalSize();
742
743 Matrix jacobian2(param2->GlobalSize(), param2->LocalSize());
744 EXPECT_TRUE(param2->ComputeJacobian(&x[x_cursor], jacobian2.data()));
745 jacobian.block(x_cursor, delta_cursor,
746 param2->GlobalSize(),
747 param2->LocalSize())
748 -= jacobian2;
749 x_cursor += param2->GlobalSize();
750 delta_cursor += param2->LocalSize();
751
752 Matrix jacobian3(param3->GlobalSize(), param3->LocalSize());
753 EXPECT_TRUE(param3->ComputeJacobian(&x[x_cursor], jacobian3.data()));
754 jacobian.block(x_cursor, delta_cursor,
755 param3->GlobalSize(),
756 param3->LocalSize())
757 -= jacobian3;
758 x_cursor += param3->GlobalSize();
759 delta_cursor += param3->LocalSize();
760
761 Matrix jacobian4(param4->GlobalSize(), param4->LocalSize());
762 EXPECT_TRUE(param4->ComputeJacobian(&x[x_cursor], jacobian4.data()));
763 jacobian.block(x_cursor, delta_cursor,
764 param4->GlobalSize(),
765 param4->LocalSize())
766 -= jacobian4;
767 x_cursor += param4->GlobalSize();
768 delta_cursor += param4->LocalSize();
769
770 EXPECT_NEAR(jacobian.norm(), 0.0, std::numeric_limits<double>::epsilon());
771}
772
773} // namespace internal
774} // namespace ceres