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