blob: 0950c0e24e2fd158f899c7772c553c90c02d9776 [file] [log] [blame]
Sameer Agarwalfa015192012-06-11 14:21:42 -07001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2012 Google Inc. All rights reserved.
3// http://code.google.com/p/ceres-solver/
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/dogleg_strategy.h"
32
33#include <cmath>
34#include "Eigen/Core"
35#include <glog/logging.h>
36#include "ceres/array_utils.h"
37#include "ceres/internal/eigen.h"
38#include "ceres/linear_solver.h"
39#include "ceres/sparse_matrix.h"
40#include "ceres/trust_region_strategy.h"
41#include "ceres/types.h"
42
43namespace ceres {
44namespace internal {
45namespace {
46const double kMaxMu = 1.0;
47const double kMinMu = 1e-8;
48}
49
50DoglegStrategy::DoglegStrategy(const TrustRegionStrategy::Options& options)
51 : linear_solver_(options.linear_solver),
52 radius_(options.initial_radius),
53 max_radius_(options.max_radius),
54 min_diagonal_(options.lm_min_diagonal),
55 max_diagonal_(options.lm_max_diagonal),
56 mu_(kMinMu),
57 min_mu_(kMinMu),
58 max_mu_(kMaxMu),
59 mu_increase_factor_(10.0),
60 increase_threshold_(0.75),
61 decrease_threshold_(0.25),
62 dogleg_step_norm_(0.0),
63 reuse_(false) {
64 CHECK_NOTNULL(linear_solver_);
65 CHECK_GT(min_diagonal_, 0.0);
66 CHECK_LT(min_diagonal_, max_diagonal_);
67 CHECK_GT(max_radius_, 0.0);
68}
69
70// If the reuse_ flag is not set, then the Cauchy point (scaled
71// gradient) and the new Gauss-Newton step are computed from
72// scratch. The Dogleg step is then computed as interpolation of these
73// two vectors.
74LinearSolver::Summary DoglegStrategy::ComputeStep(
75 const TrustRegionStrategy::PerSolveOptions& per_solve_options,
76 SparseMatrix* jacobian,
77 const double* residuals,
78 double* step) {
79 CHECK_NOTNULL(jacobian);
80 CHECK_NOTNULL(residuals);
81 CHECK_NOTNULL(step);
82
83 const int n = jacobian->num_cols();
84 if (reuse_) {
85 // Gauss-Newton and gradient vectors are always available, only a
86 // new interpolant need to be computed.
87 ComputeDoglegStep(step);
88 LinearSolver::Summary linear_solver_summary;
89 linear_solver_summary.num_iterations = 0;
90 linear_solver_summary.termination_type = TOLERANCE;
91 return linear_solver_summary;
92 }
93
94 reuse_ = true;
95 // Check that we have the storage needed to hold the various
96 // temporary vectors.
97 if (diagonal_.rows() != n) {
98 diagonal_.resize(n, 1);
99 gradient_.resize(n, 1);
100 gauss_newton_step_.resize(n, 1);
101 }
102
103 // Vector used to form the diagonal matrix that is used to
104 // regularize the Gauss-Newton solve.
105 jacobian->SquaredColumnNorm(diagonal_.data());
106 for (int i = 0; i < n; ++i) {
107 diagonal_[i] = min(max(diagonal_[i], min_diagonal_), max_diagonal_);
108 }
109
110 gradient_.setZero();
111 jacobian->LeftMultiply(residuals, gradient_.data());
112
113 // alpha * gradient is the Cauchy point.
114 Vector Jg(jacobian->num_rows());
115 Jg.setZero();
116 jacobian->RightMultiply(gradient_.data(), Jg.data());
117 alpha_ = gradient_.squaredNorm() / Jg.squaredNorm();
118
119 LinearSolver::Summary linear_solver_summary =
120 ComputeGaussNewtonStep(jacobian, residuals);
121
122 // Interpolate the Cauchy point and the Gauss-Newton step.
123 if (linear_solver_summary.termination_type != FAILURE) {
124 ComputeDoglegStep(step);
125 }
126
127 return linear_solver_summary;
128}
129
130void DoglegStrategy::ComputeDoglegStep(double* dogleg) {
131 VectorRef dogleg_step(dogleg, gradient_.rows());
132
133 // Case 1. The Gauss-Newton step lies inside the trust region, and
134 // is therefore the optimal solution to the trust-region problem.
135 const double gradient_norm = gradient_.norm();
136 const double gauss_newton_norm = gauss_newton_step_.norm();
137 if (gauss_newton_norm <= radius_) {
138 dogleg_step = gauss_newton_step_;
139 dogleg_step_norm_ = gauss_newton_norm;
140 VLOG(3) << "GaussNewton step size: " << dogleg_step_norm_
141 << " radius: " << radius_;
142 return;
143 }
144
145 // Case 2. The Cauchy point and the Gauss-Newton steps lie outside
146 // the trust region. Rescale the Cauchy point to the trust region
147 // and return.
148 if (gradient_norm * alpha_ >= radius_) {
149 dogleg_step = (radius_ / gradient_norm) * gradient_;
150 dogleg_step_norm_ = radius_;
151 VLOG(3) << "Cauchy step size: " << dogleg_step_norm_
152 << " radius: " << radius_;
153 return;
154 }
155
156 // Case 3. The Cauchy point is inside the trust region and the
157 // Gauss-Newton step is outside. Compute the line joining the two
158 // points and the point on it which intersects the trust region
159 // boundary.
160
161 // a = alpha * gradient
162 // b = gauss_newton_step
163 const double b_dot_a = alpha_ * gradient_.dot(gauss_newton_step_);
164 const double a_squared_norm = pow(alpha_ * gradient_norm, 2.0);
165 const double b_minus_a_squared_norm =
166 a_squared_norm - 2 * b_dot_a + pow(gauss_newton_norm, 2);
167
168 // c = a' (b - a)
169 // = alpha * gradient' gauss_newton_step - alpha^2 |gradient|^2
170 const double c = b_dot_a - a_squared_norm;
171 const double d = sqrt(c * c + b_minus_a_squared_norm *
172 (pow(radius_, 2.0) - a_squared_norm));
173
174 double beta =
175 (c <= 0)
176 ? (d - c) / b_minus_a_squared_norm
177 : (radius_ * radius_ - a_squared_norm) / (d + c);
178 dogleg_step = (alpha_ * (1.0 - beta)) * gradient_ + beta * gauss_newton_step_;
179 dogleg_step_norm_ = dogleg_step.norm();
180 VLOG(3) << "Dogleg step size: " << dogleg_step_norm_
181 << " radius: " << radius_;
182}
183
184LinearSolver::Summary DoglegStrategy::ComputeGaussNewtonStep(
185 SparseMatrix* jacobian,
186 const double* residuals) {
187 const int n = jacobian->num_cols();
188 LinearSolver::Summary linear_solver_summary;
189 linear_solver_summary.termination_type = FAILURE;
190
191 // The Jacobian matrix is often quite poorly conditioned. Thus it is
192 // necessary to add a diagonal matrix at the bottom to prevent the
193 // linear solver from failing.
194 //
195 // We do this by computing the same diagonal matrix as the one used
196 // by Levenberg-Marquardt (other choices are possible), and scaling
197 // it by a small constant (independent of the trust region radius).
198 //
199 // If the solve fails, the multiplier to the diagonal is increased
200 // up to max_mu_ by a factor of mu_increase_factor_ every time. If
201 // the linear solver is still not successful, the strategy returns
202 // with FAILURE.
203 //
204 // Next time when a new Gauss-Newton step is requested, the
205 // multiplier starts out from the last successful solve.
206 //
207 // When a step is declared successful, the multiplier is decreased
208 // by half of mu_increase_factor_.
209 while (mu_ < max_mu_) {
210 // Dogleg, as far as I (sameeragarwal) understand it, requires a
211 // reasonably good estimate of the Gauss-Newton step. This means
212 // that we need to solve the normal equations more or less
213 // exactly. This is reflected in the values of the tolerances set
214 // below.
215 //
216 // For now, this strategy should only be used with exact
217 // factorization based solvers, for which these tolerances are
218 // automatically satisfied.
219 //
220 // The right way to combine inexact solves with trust region
221 // methods is to use Stiehaug's method.
222 LinearSolver::PerSolveOptions solve_options;
223 solve_options.q_tolerance = 0.0;
224 solve_options.r_tolerance = 0.0;
225
226 lm_diagonal_ = (diagonal_ * mu_).array().sqrt();
227 solve_options.D = lm_diagonal_.data();
228
229 InvalidateArray(n, gauss_newton_step_.data());
230 linear_solver_summary = linear_solver_->Solve(jacobian,
231 residuals,
232 solve_options,
233 gauss_newton_step_.data());
234
235 if (linear_solver_summary.termination_type == FAILURE ||
236 !IsArrayValid(n, gauss_newton_step_.data())) {
237 mu_ *= mu_increase_factor_;
238 VLOG(2) << "Increasing mu " << mu_;
239 linear_solver_summary.termination_type = FAILURE;
240 continue;
241 }
242 break;
243 }
244
245 return linear_solver_summary;
246}
247
248void DoglegStrategy::StepAccepted(double step_quality) {
249 CHECK_GT(step_quality, 0.0);
250 if (step_quality < decrease_threshold_) {
251 radius_ *= 0.5;
252 return;
253 }
254
255 if (step_quality > increase_threshold_) {
256 radius_ = max(radius_, 3.0 * dogleg_step_norm_);
257 }
258
259 // Reduce the regularization multiplier, in the hope that whatever
260 // was causing the rank deficiency has gone away and we can return
261 // to doing a pure Gauss-Newton solve.
262 mu_ = max(min_mu_, 2.0 * mu_ / mu_increase_factor_ );
263 reuse_ = false;
264}
265
266void DoglegStrategy::StepRejected(double step_quality) {
267 radius_ *= 0.5;
268 reuse_ = true;
269}
270
271void DoglegStrategy::StepIsInvalid() {
272 mu_ *= mu_increase_factor_;
273 reuse_ = false;
274}
275
276double DoglegStrategy::Radius() const {
277 return radius_;
278}
279
280} // namespace internal
281} // namespace ceres