blob: a7fe4a1afe5c418cdeee500fe00d9b64c8635731 [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 "ceres/local_parameterization.h"
32
33#include <algorithm>
34#include "Eigen/Geometry"
35#include "ceres/householder_vector.h"
36#include "ceres/internal/eigen.h"
37#include "ceres/internal/fixed_array.h"
38#include "ceres/rotation.h"
39#include "glog/logging.h"
40
41namespace ceres {
42
43using std::vector;
44
45LocalParameterization::~LocalParameterization() {
46}
47
48bool LocalParameterization::MultiplyByJacobian(const double* x,
49 const int num_rows,
50 const double* global_matrix,
51 double* local_matrix) const {
52 if (LocalSize() == 0) {
53 return true;
54 }
55
56 Matrix jacobian(GlobalSize(), LocalSize());
57 if (!ComputeJacobian(x, jacobian.data())) {
58 return false;
59 }
60
61 MatrixRef(local_matrix, num_rows, LocalSize()) =
62 ConstMatrixRef(global_matrix, num_rows, GlobalSize()) * jacobian;
63 return true;
64}
65
66IdentityParameterization::IdentityParameterization(const int size)
67 : size_(size) {
68 CHECK_GT(size, 0);
69}
70
71bool IdentityParameterization::Plus(const double* x,
72 const double* delta,
73 double* x_plus_delta) const {
74 VectorRef(x_plus_delta, size_) =
75 ConstVectorRef(x, size_) + ConstVectorRef(delta, size_);
76 return true;
77}
78
79bool IdentityParameterization::ComputeJacobian(const double* x,
80 double* jacobian) const {
81 MatrixRef(jacobian, size_, size_).setIdentity();
82 return true;
83}
84
85bool IdentityParameterization::MultiplyByJacobian(const double* x,
86 const int num_cols,
87 const double* global_matrix,
88 double* local_matrix) const {
89 std::copy(global_matrix,
90 global_matrix + num_cols * GlobalSize(),
91 local_matrix);
92 return true;
93}
94
95SubsetParameterization::SubsetParameterization(
96 int size, const vector<int>& constant_parameters)
97 : local_size_(size - constant_parameters.size()), constancy_mask_(size, 0) {
98 vector<int> constant = constant_parameters;
99 std::sort(constant.begin(), constant.end());
100 CHECK_GE(constant.front(), 0) << "Indices indicating constant parameter must "
101 "be greater than equal to zero.";
102 CHECK_LT(constant.back(), size)
103 << "Indices indicating constant parameter must be less than the size "
104 << "of the parameter block.";
105 CHECK(std::adjacent_find(constant.begin(), constant.end()) == constant.end())
106 << "The set of constant parameters cannot contain duplicates";
107 for (int i = 0; i < constant_parameters.size(); ++i) {
108 constancy_mask_[constant_parameters[i]] = 1;
109 }
110}
111
112bool SubsetParameterization::Plus(const double* x,
113 const double* delta,
114 double* x_plus_delta) const {
115 const int global_size = GlobalSize();
116 for (int i = 0, j = 0; i < global_size; ++i) {
117 if (constancy_mask_[i]) {
118 x_plus_delta[i] = x[i];
119 } else {
120 x_plus_delta[i] = x[i] + delta[j++];
121 }
122 }
123 return true;
124}
125
126bool SubsetParameterization::ComputeJacobian(const double* x,
127 double* jacobian) const {
128 if (local_size_ == 0) {
129 return true;
130 }
131
132 const int global_size = GlobalSize();
133 MatrixRef m(jacobian, global_size, local_size_);
134 m.setZero();
135 for (int i = 0, j = 0; i < global_size; ++i) {
136 if (!constancy_mask_[i]) {
137 m(i, j++) = 1.0;
138 }
139 }
140 return true;
141}
142
143bool SubsetParameterization::MultiplyByJacobian(const double* x,
144 const int num_rows,
145 const double* global_matrix,
146 double* local_matrix) const {
147 if (local_size_ == 0) {
148 return true;
149 }
150
151 const int global_size = GlobalSize();
152 for (int row = 0; row < num_rows; ++row) {
153 for (int col = 0, j = 0; col < global_size; ++col) {
154 if (!constancy_mask_[col]) {
155 local_matrix[row * local_size_ + j++] =
156 global_matrix[row * global_size + col];
157 }
158 }
159 }
160 return true;
161}
162
163bool QuaternionParameterization::Plus(const double* x,
164 const double* delta,
165 double* x_plus_delta) const {
166 const double norm_delta =
167 sqrt(delta[0] * delta[0] + delta[1] * delta[1] + delta[2] * delta[2]);
168 if (norm_delta > 0.0) {
169 const double sin_delta_by_delta = (sin(norm_delta) / norm_delta);
170 double q_delta[4];
171 q_delta[0] = cos(norm_delta);
172 q_delta[1] = sin_delta_by_delta * delta[0];
173 q_delta[2] = sin_delta_by_delta * delta[1];
174 q_delta[3] = sin_delta_by_delta * delta[2];
175 QuaternionProduct(q_delta, x, x_plus_delta);
176 } else {
177 for (int i = 0; i < 4; ++i) {
178 x_plus_delta[i] = x[i];
179 }
180 }
181 return true;
182}
183
184bool QuaternionParameterization::ComputeJacobian(const double* x,
185 double* jacobian) const {
186 jacobian[0] = -x[1]; jacobian[1] = -x[2]; jacobian[2] = -x[3]; // NOLINT
187 jacobian[3] = x[0]; jacobian[4] = x[3]; jacobian[5] = -x[2]; // NOLINT
188 jacobian[6] = -x[3]; jacobian[7] = x[0]; jacobian[8] = x[1]; // NOLINT
189 jacobian[9] = x[2]; jacobian[10] = -x[1]; jacobian[11] = x[0]; // NOLINT
190 return true;
191}
192
193bool EigenQuaternionParameterization::Plus(const double* x_ptr,
194 const double* delta,
195 double* x_plus_delta_ptr) const {
196 Eigen::Map<Eigen::Quaterniond> x_plus_delta(x_plus_delta_ptr);
197 Eigen::Map<const Eigen::Quaterniond> x(x_ptr);
198
199 const double norm_delta =
200 sqrt(delta[0] * delta[0] + delta[1] * delta[1] + delta[2] * delta[2]);
201 if (norm_delta > 0.0) {
202 const double sin_delta_by_delta = sin(norm_delta) / norm_delta;
203
204 // Note, in the constructor w is first.
205 Eigen::Quaterniond delta_q(cos(norm_delta),
206 sin_delta_by_delta * delta[0],
207 sin_delta_by_delta * delta[1],
208 sin_delta_by_delta * delta[2]);
209 x_plus_delta = delta_q * x;
210 } else {
211 x_plus_delta = x;
212 }
213
214 return true;
215}
216
217bool EigenQuaternionParameterization::ComputeJacobian(const double* x,
218 double* jacobian) const {
219 jacobian[0] = x[3]; jacobian[1] = x[2]; jacobian[2] = -x[1]; // NOLINT
220 jacobian[3] = -x[2]; jacobian[4] = x[3]; jacobian[5] = x[0]; // NOLINT
221 jacobian[6] = x[1]; jacobian[7] = -x[0]; jacobian[8] = x[3]; // NOLINT
222 jacobian[9] = -x[0]; jacobian[10] = -x[1]; jacobian[11] = -x[2]; // NOLINT
223 return true;
224}
225
226HomogeneousVectorParameterization::HomogeneousVectorParameterization(int size)
227 : size_(size) {
228 CHECK_GT(size_, 1) << "The size of the homogeneous vector needs to be "
229 << "greater than 1.";
230}
231
232bool HomogeneousVectorParameterization::Plus(const double* x_ptr,
233 const double* delta_ptr,
234 double* x_plus_delta_ptr) const {
235 ConstVectorRef x(x_ptr, size_);
236 ConstVectorRef delta(delta_ptr, size_ - 1);
237 VectorRef x_plus_delta(x_plus_delta_ptr, size_);
238
239 const double norm_delta = delta.norm();
240
241 if (norm_delta == 0.0) {
242 x_plus_delta = x;
243 return true;
244 }
245
246 // Map the delta from the minimum representation to the over parameterized
247 // homogeneous vector. See section A6.9.2 on page 624 of Hartley & Zisserman
248 // (2nd Edition) for a detailed description. Note there is a typo on Page
249 // 625, line 4 so check the book errata.
250 const double norm_delta_div_2 = 0.5 * norm_delta;
251 const double sin_delta_by_delta = sin(norm_delta_div_2) /
252 norm_delta_div_2;
253
254 Vector y(size_);
255 y.head(size_ - 1) = 0.5 * sin_delta_by_delta * delta;
256 y(size_ - 1) = cos(norm_delta_div_2);
257
258 Vector v(size_);
259 double beta;
260 internal::ComputeHouseholderVector<double>(x, &v, &beta);
261
262 // Apply the delta update to remain on the unit sphere. See section A6.9.3
263 // on page 625 of Hartley & Zisserman (2nd Edition) for a detailed
264 // description.
265 x_plus_delta = x.norm() * (y - v * (beta * (v.transpose() * y)));
266
267 return true;
268}
269
270bool HomogeneousVectorParameterization::ComputeJacobian(
271 const double* x_ptr, double* jacobian_ptr) const {
272 ConstVectorRef x(x_ptr, size_);
273 MatrixRef jacobian(jacobian_ptr, size_, size_ - 1);
274
275 Vector v(size_);
276 double beta;
277 internal::ComputeHouseholderVector<double>(x, &v, &beta);
278
279 // The Jacobian is equal to J = 0.5 * H.leftCols(size_ - 1) where H is the
280 // Householder matrix (H = I - beta * v * v').
281 for (int i = 0; i < size_ - 1; ++i) {
282 jacobian.col(i) = -0.5 * beta * v(i) * v;
283 jacobian.col(i)(i) += 0.5;
284 }
285 jacobian *= x.norm();
286
287 return true;
288}
289
290ProductParameterization::ProductParameterization(
291 LocalParameterization* local_param1,
292 LocalParameterization* local_param2) {
293 local_params_.push_back(local_param1);
294 local_params_.push_back(local_param2);
295 Init();
296}
297
298ProductParameterization::ProductParameterization(
299 LocalParameterization* local_param1,
300 LocalParameterization* local_param2,
301 LocalParameterization* local_param3) {
302 local_params_.push_back(local_param1);
303 local_params_.push_back(local_param2);
304 local_params_.push_back(local_param3);
305 Init();
306}
307
308ProductParameterization::ProductParameterization(
309 LocalParameterization* local_param1,
310 LocalParameterization* local_param2,
311 LocalParameterization* local_param3,
312 LocalParameterization* local_param4) {
313 local_params_.push_back(local_param1);
314 local_params_.push_back(local_param2);
315 local_params_.push_back(local_param3);
316 local_params_.push_back(local_param4);
317 Init();
318}
319
320ProductParameterization::~ProductParameterization() {
321 for (int i = 0; i < local_params_.size(); ++i) {
322 delete local_params_[i];
323 }
324}
325
326void ProductParameterization::Init() {
327 global_size_ = 0;
328 local_size_ = 0;
329 buffer_size_ = 0;
330 for (int i = 0; i < local_params_.size(); ++i) {
331 const LocalParameterization* param = local_params_[i];
332 buffer_size_ = std::max(buffer_size_,
333 param->LocalSize() * param->GlobalSize());
334 global_size_ += param->GlobalSize();
335 local_size_ += param->LocalSize();
336 }
337}
338
339bool ProductParameterization::Plus(const double* x,
340 const double* delta,
341 double* x_plus_delta) const {
342 int x_cursor = 0;
343 int delta_cursor = 0;
344 for (int i = 0; i < local_params_.size(); ++i) {
345 const LocalParameterization* param = local_params_[i];
346 if (!param->Plus(x + x_cursor,
347 delta + delta_cursor,
348 x_plus_delta + x_cursor)) {
349 return false;
350 }
351 delta_cursor += param->LocalSize();
352 x_cursor += param->GlobalSize();
353 }
354
355 return true;
356}
357
358bool ProductParameterization::ComputeJacobian(const double* x,
359 double* jacobian_ptr) const {
360 MatrixRef jacobian(jacobian_ptr, GlobalSize(), LocalSize());
361 jacobian.setZero();
362 internal::FixedArray<double> buffer(buffer_size_);
363
364 int x_cursor = 0;
365 int delta_cursor = 0;
366 for (int i = 0; i < local_params_.size(); ++i) {
367 const LocalParameterization* param = local_params_[i];
368 const int local_size = param->LocalSize();
369 const int global_size = param->GlobalSize();
370
371 if (!param->ComputeJacobian(x + x_cursor, buffer.get())) {
372 return false;
373 }
374 jacobian.block(x_cursor, delta_cursor, global_size, local_size)
375 = MatrixRef(buffer.get(), global_size, local_size);
376
377 delta_cursor += local_size;
378 x_cursor += global_size;
379 }
380
381 return true;
382}
383
384} // namespace ceres