Imported Upstream version ceres 1.13.0
[platform/upstream/ceres-solver.git] / internal / ceres / levenberg_marquardt_strategy.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/levenberg_marquardt_strategy.h"
32
33 #include <cmath>
34 #include "Eigen/Core"
35 #include "ceres/array_utils.h"
36 #include "ceres/internal/eigen.h"
37 #include "ceres/linear_least_squares_problems.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 #include "glog/logging.h"
43
44 namespace ceres {
45 namespace internal {
46
47 LevenbergMarquardtStrategy::LevenbergMarquardtStrategy(
48     const TrustRegionStrategy::Options& options)
49     : linear_solver_(options.linear_solver),
50       radius_(options.initial_radius),
51       max_radius_(options.max_radius),
52       min_diagonal_(options.min_lm_diagonal),
53       max_diagonal_(options.max_lm_diagonal),
54       decrease_factor_(2.0),
55       reuse_diagonal_(false) {
56   CHECK_NOTNULL(linear_solver_);
57   CHECK_GT(min_diagonal_, 0.0);
58   CHECK_LE(min_diagonal_, max_diagonal_);
59   CHECK_GT(max_radius_, 0.0);
60 }
61
62 LevenbergMarquardtStrategy::~LevenbergMarquardtStrategy() {
63 }
64
65 TrustRegionStrategy::Summary LevenbergMarquardtStrategy::ComputeStep(
66     const TrustRegionStrategy::PerSolveOptions& per_solve_options,
67     SparseMatrix* jacobian,
68     const double* residuals,
69     double* step) {
70   CHECK_NOTNULL(jacobian);
71   CHECK_NOTNULL(residuals);
72   CHECK_NOTNULL(step);
73
74   const int num_parameters = jacobian->num_cols();
75   if (!reuse_diagonal_) {
76     if (diagonal_.rows() != num_parameters) {
77       diagonal_.resize(num_parameters, 1);
78     }
79
80     jacobian->SquaredColumnNorm(diagonal_.data());
81     for (int i = 0; i < num_parameters; ++i) {
82       diagonal_[i] = std::min(std::max(diagonal_[i], min_diagonal_),
83                               max_diagonal_);
84     }
85   }
86
87   lm_diagonal_ = (diagonal_ / radius_).array().sqrt();
88
89   LinearSolver::PerSolveOptions solve_options;
90   solve_options.D = lm_diagonal_.data();
91   solve_options.q_tolerance = per_solve_options.eta;
92   // Disable r_tolerance checking. Since we only care about
93   // termination via the q_tolerance. As Nash and Sofer show,
94   // r_tolerance based termination is essentially useless in
95   // Truncated Newton methods.
96   solve_options.r_tolerance = -1.0;
97
98   // Invalidate the output array lm_step, so that we can detect if
99   // the linear solver generated numerical garbage.  This is known
100   // to happen for the DENSE_QR and then DENSE_SCHUR solver when
101   // the Jacobin is severly rank deficient and mu is too small.
102   InvalidateArray(num_parameters, step);
103
104   // Instead of solving Jx = -r, solve Jy = r.
105   // Then x can be found as x = -y, but the inputs jacobian and residuals
106   // do not need to be modified.
107   LinearSolver::Summary linear_solver_summary =
108       linear_solver_->Solve(jacobian, residuals, solve_options, step);
109
110   if (linear_solver_summary.termination_type == LINEAR_SOLVER_FATAL_ERROR) {
111     LOG(WARNING) << "Linear solver fatal error: "
112                  << linear_solver_summary.message;
113   } else if (linear_solver_summary.termination_type == LINEAR_SOLVER_FAILURE)  {
114     LOG(WARNING) << "Linear solver failure. Failed to compute a step: "
115                  << linear_solver_summary.message;
116   } else if (!IsArrayValid(num_parameters, step)) {
117     LOG(WARNING) << "Linear solver failure. Failed to compute a finite step.";
118     linear_solver_summary.termination_type = LINEAR_SOLVER_FAILURE;
119   } else {
120     VectorRef(step, num_parameters) *= -1.0;
121   }
122   reuse_diagonal_ = true;
123
124   if (per_solve_options.dump_format_type == CONSOLE ||
125       (per_solve_options.dump_format_type != CONSOLE &&
126        !per_solve_options.dump_filename_base.empty())) {
127     if (!DumpLinearLeastSquaresProblem(per_solve_options.dump_filename_base,
128                                        per_solve_options.dump_format_type,
129                                        jacobian,
130                                        solve_options.D,
131                                        residuals,
132                                        step,
133                                        0)) {
134       LOG(ERROR) << "Unable to dump trust region problem."
135                  << " Filename base: " << per_solve_options.dump_filename_base;
136     }
137   }
138
139
140   TrustRegionStrategy::Summary summary;
141   summary.residual_norm = linear_solver_summary.residual_norm;
142   summary.num_iterations = linear_solver_summary.num_iterations;
143   summary.termination_type = linear_solver_summary.termination_type;
144   return summary;
145 }
146
147 void LevenbergMarquardtStrategy::StepAccepted(double step_quality) {
148   CHECK_GT(step_quality, 0.0);
149   radius_ = radius_ / std::max(1.0 / 3.0,
150                                1.0 - pow(2.0 * step_quality - 1.0, 3));
151   radius_ = std::min(max_radius_, radius_);
152   decrease_factor_ = 2.0;
153   reuse_diagonal_ = false;
154 }
155
156 void LevenbergMarquardtStrategy::StepRejected(double step_quality) {
157   radius_ = radius_ / decrease_factor_;
158   decrease_factor_ *= 2.0;
159   reuse_diagonal_ = true;
160 }
161
162 double LevenbergMarquardtStrategy::Radius() const {
163   return radius_;
164 }
165
166 }  // namespace internal
167 }  // namespace ceres