blob: 62947f06fccdea3d45edc72c48711f6a83f7d67d [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>
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080034
Austin Schuh70cc9552019-01-21 19:46:48 -080035#include "Eigen/Geometry"
Austin Schuh70cc9552019-01-21 19:46:48 -080036#include "ceres/internal/eigen.h"
37#include "ceres/internal/fixed_array.h"
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080038#include "ceres/internal/householder_vector.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080039#include "ceres/rotation.h"
40#include "glog/logging.h"
41
42namespace ceres {
43
44using std::vector;
45
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080046LocalParameterization::~LocalParameterization() {}
Austin Schuh70cc9552019-01-21 19:46:48 -080047
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 {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080089 std::copy(
90 global_matrix, global_matrix + num_cols * GlobalSize(), local_matrix);
Austin Schuh70cc9552019-01-21 19:46:48 -080091 return true;
92}
93
94SubsetParameterization::SubsetParameterization(
95 int size, const vector<int>& constant_parameters)
96 : local_size_(size - constant_parameters.size()), constancy_mask_(size, 0) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080097 if (constant_parameters.empty()) {
98 return;
99 }
100
Austin Schuh70cc9552019-01-21 19:46:48 -0800101 vector<int> constant = constant_parameters;
102 std::sort(constant.begin(), constant.end());
103 CHECK_GE(constant.front(), 0) << "Indices indicating constant parameter must "
104 "be greater than equal to zero.";
105 CHECK_LT(constant.back(), size)
106 << "Indices indicating constant parameter must be less than the size "
107 << "of the parameter block.";
108 CHECK(std::adjacent_find(constant.begin(), constant.end()) == constant.end())
109 << "The set of constant parameters cannot contain duplicates";
110 for (int i = 0; i < constant_parameters.size(); ++i) {
111 constancy_mask_[constant_parameters[i]] = 1;
112 }
113}
114
115bool SubsetParameterization::Plus(const double* x,
116 const double* delta,
117 double* x_plus_delta) const {
118 const int global_size = GlobalSize();
119 for (int i = 0, j = 0; i < global_size; ++i) {
120 if (constancy_mask_[i]) {
121 x_plus_delta[i] = x[i];
122 } else {
123 x_plus_delta[i] = x[i] + delta[j++];
124 }
125 }
126 return true;
127}
128
129bool SubsetParameterization::ComputeJacobian(const double* x,
130 double* jacobian) const {
131 if (local_size_ == 0) {
132 return true;
133 }
134
135 const int global_size = GlobalSize();
136 MatrixRef m(jacobian, global_size, local_size_);
137 m.setZero();
138 for (int i = 0, j = 0; i < global_size; ++i) {
139 if (!constancy_mask_[i]) {
140 m(i, j++) = 1.0;
141 }
142 }
143 return true;
144}
145
146bool SubsetParameterization::MultiplyByJacobian(const double* x,
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800147 const int num_cols,
148 const double* global_matrix,
149 double* local_matrix) const {
Austin Schuh70cc9552019-01-21 19:46:48 -0800150 if (local_size_ == 0) {
151 return true;
152 }
153
154 const int global_size = GlobalSize();
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800155 for (int col = 0; col < num_cols; ++col) {
156 for (int i = 0, j = 0; i < global_size; ++i) {
157 if (!constancy_mask_[i]) {
158 local_matrix[col * local_size_ + j++] =
159 global_matrix[col * global_size + i];
Austin Schuh70cc9552019-01-21 19:46:48 -0800160 }
161 }
162 }
163 return true;
164}
165
166bool QuaternionParameterization::Plus(const double* x,
167 const double* delta,
168 double* x_plus_delta) const {
169 const double norm_delta =
170 sqrt(delta[0] * delta[0] + delta[1] * delta[1] + delta[2] * delta[2]);
171 if (norm_delta > 0.0) {
172 const double sin_delta_by_delta = (sin(norm_delta) / norm_delta);
173 double q_delta[4];
174 q_delta[0] = cos(norm_delta);
175 q_delta[1] = sin_delta_by_delta * delta[0];
176 q_delta[2] = sin_delta_by_delta * delta[1];
177 q_delta[3] = sin_delta_by_delta * delta[2];
178 QuaternionProduct(q_delta, x, x_plus_delta);
179 } else {
180 for (int i = 0; i < 4; ++i) {
181 x_plus_delta[i] = x[i];
182 }
183 }
184 return true;
185}
186
187bool QuaternionParameterization::ComputeJacobian(const double* x,
188 double* jacobian) const {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800189 // clang-format off
190 jacobian[0] = -x[1]; jacobian[1] = -x[2]; jacobian[2] = -x[3];
191 jacobian[3] = x[0]; jacobian[4] = x[3]; jacobian[5] = -x[2];
192 jacobian[6] = -x[3]; jacobian[7] = x[0]; jacobian[8] = x[1];
193 jacobian[9] = x[2]; jacobian[10] = -x[1]; jacobian[11] = x[0];
194 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800195 return true;
196}
197
198bool EigenQuaternionParameterization::Plus(const double* x_ptr,
199 const double* delta,
200 double* x_plus_delta_ptr) const {
201 Eigen::Map<Eigen::Quaterniond> x_plus_delta(x_plus_delta_ptr);
202 Eigen::Map<const Eigen::Quaterniond> x(x_ptr);
203
204 const double norm_delta =
205 sqrt(delta[0] * delta[0] + delta[1] * delta[1] + delta[2] * delta[2]);
206 if (norm_delta > 0.0) {
207 const double sin_delta_by_delta = sin(norm_delta) / norm_delta;
208
209 // Note, in the constructor w is first.
210 Eigen::Quaterniond delta_q(cos(norm_delta),
211 sin_delta_by_delta * delta[0],
212 sin_delta_by_delta * delta[1],
213 sin_delta_by_delta * delta[2]);
214 x_plus_delta = delta_q * x;
215 } else {
216 x_plus_delta = x;
217 }
218
219 return true;
220}
221
222bool EigenQuaternionParameterization::ComputeJacobian(const double* x,
223 double* jacobian) const {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800224 // clang-format off
225 jacobian[0] = x[3]; jacobian[1] = x[2]; jacobian[2] = -x[1];
226 jacobian[3] = -x[2]; jacobian[4] = x[3]; jacobian[5] = x[0];
227 jacobian[6] = x[1]; jacobian[7] = -x[0]; jacobian[8] = x[3];
228 jacobian[9] = -x[0]; jacobian[10] = -x[1]; jacobian[11] = -x[2];
229 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800230 return true;
231}
232
233HomogeneousVectorParameterization::HomogeneousVectorParameterization(int size)
234 : size_(size) {
235 CHECK_GT(size_, 1) << "The size of the homogeneous vector needs to be "
236 << "greater than 1.";
237}
238
239bool HomogeneousVectorParameterization::Plus(const double* x_ptr,
240 const double* delta_ptr,
241 double* x_plus_delta_ptr) const {
242 ConstVectorRef x(x_ptr, size_);
243 ConstVectorRef delta(delta_ptr, size_ - 1);
244 VectorRef x_plus_delta(x_plus_delta_ptr, size_);
245
246 const double norm_delta = delta.norm();
247
248 if (norm_delta == 0.0) {
249 x_plus_delta = x;
250 return true;
251 }
252
253 // Map the delta from the minimum representation to the over parameterized
254 // homogeneous vector. See section A6.9.2 on page 624 of Hartley & Zisserman
255 // (2nd Edition) for a detailed description. Note there is a typo on Page
256 // 625, line 4 so check the book errata.
257 const double norm_delta_div_2 = 0.5 * norm_delta;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800258 const double sin_delta_by_delta =
259 std::sin(norm_delta_div_2) / norm_delta_div_2;
Austin Schuh70cc9552019-01-21 19:46:48 -0800260
261 Vector y(size_);
262 y.head(size_ - 1) = 0.5 * sin_delta_by_delta * delta;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800263 y(size_ - 1) = std::cos(norm_delta_div_2);
Austin Schuh70cc9552019-01-21 19:46:48 -0800264
265 Vector v(size_);
266 double beta;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800267
268 // NOTE: The explicit template arguments are needed here because
269 // ComputeHouseholderVector is templated and some versions of MSVC
270 // have trouble deducing the type of v automatically.
271 internal::ComputeHouseholderVector<ConstVectorRef, double, Eigen::Dynamic>(
272 x, &v, &beta);
Austin Schuh70cc9552019-01-21 19:46:48 -0800273
274 // Apply the delta update to remain on the unit sphere. See section A6.9.3
275 // on page 625 of Hartley & Zisserman (2nd Edition) for a detailed
276 // description.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800277 x_plus_delta = x.norm() * (y - v * (beta * (v.transpose() * y)));
Austin Schuh70cc9552019-01-21 19:46:48 -0800278
279 return true;
280}
281
282bool HomogeneousVectorParameterization::ComputeJacobian(
283 const double* x_ptr, double* jacobian_ptr) const {
284 ConstVectorRef x(x_ptr, size_);
285 MatrixRef jacobian(jacobian_ptr, size_, size_ - 1);
286
287 Vector v(size_);
288 double beta;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800289
290 // NOTE: The explicit template arguments are needed here because
291 // ComputeHouseholderVector is templated and some versions of MSVC
292 // have trouble deducing the type of v automatically.
293 internal::ComputeHouseholderVector<ConstVectorRef, double, Eigen::Dynamic>(
294 x, &v, &beta);
Austin Schuh70cc9552019-01-21 19:46:48 -0800295
296 // The Jacobian is equal to J = 0.5 * H.leftCols(size_ - 1) where H is the
297 // Householder matrix (H = I - beta * v * v').
298 for (int i = 0; i < size_ - 1; ++i) {
299 jacobian.col(i) = -0.5 * beta * v(i) * v;
300 jacobian.col(i)(i) += 0.5;
301 }
302 jacobian *= x.norm();
303
304 return true;
305}
306
Austin Schuh70cc9552019-01-21 19:46:48 -0800307bool ProductParameterization::Plus(const double* x,
308 const double* delta,
309 double* x_plus_delta) const {
310 int x_cursor = 0;
311 int delta_cursor = 0;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800312 for (const auto& param : local_params_) {
313 if (!param->Plus(
314 x + x_cursor, delta + delta_cursor, x_plus_delta + x_cursor)) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800315 return false;
316 }
317 delta_cursor += param->LocalSize();
318 x_cursor += param->GlobalSize();
319 }
320
321 return true;
322}
323
324bool ProductParameterization::ComputeJacobian(const double* x,
325 double* jacobian_ptr) const {
326 MatrixRef jacobian(jacobian_ptr, GlobalSize(), LocalSize());
327 jacobian.setZero();
328 internal::FixedArray<double> buffer(buffer_size_);
329
330 int x_cursor = 0;
331 int delta_cursor = 0;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800332 for (const auto& param : local_params_) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800333 const int local_size = param->LocalSize();
334 const int global_size = param->GlobalSize();
335
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800336 if (!param->ComputeJacobian(x + x_cursor, buffer.data())) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800337 return false;
338 }
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800339 jacobian.block(x_cursor, delta_cursor, global_size, local_size) =
340 MatrixRef(buffer.data(), global_size, local_size);
Austin Schuh70cc9552019-01-21 19:46:48 -0800341
342 delta_cursor += local_size;
343 x_cursor += global_size;
344 }
345
346 return true;
347}
348
349} // namespace ceres