blob: bcf05b3ddfb26a60b9559d4a633e79b0da719b1d [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2016 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/trust_region_minimizer.h"
32
33#include <algorithm>
34#include <cmath>
35#include <cstdlib>
36#include <cstring>
Austin Schuh70cc9552019-01-21 19:46:48 -080037#include <limits>
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080038#include <memory>
Austin Schuh70cc9552019-01-21 19:46:48 -080039#include <string>
40#include <vector>
41
42#include "Eigen/Core"
43#include "ceres/array_utils.h"
44#include "ceres/coordinate_descent_minimizer.h"
45#include "ceres/evaluator.h"
46#include "ceres/file.h"
47#include "ceres/line_search.h"
48#include "ceres/stringprintf.h"
49#include "ceres/types.h"
50#include "ceres/wall_time.h"
51#include "glog/logging.h"
52
53// Helper macro to simplify some of the control flow.
54#define RETURN_IF_ERROR_AND_LOG(expr) \
55 do { \
56 if (!(expr)) { \
57 LOG(ERROR) << "Terminating: " << solver_summary_->message; \
58 return; \
59 } \
60 } while (0)
61
62namespace ceres {
63namespace internal {
64
65TrustRegionMinimizer::~TrustRegionMinimizer() {}
66
67void TrustRegionMinimizer::Minimize(const Minimizer::Options& options,
68 double* parameters,
69 Solver::Summary* solver_summary) {
70 start_time_in_secs_ = WallTimeInSeconds();
71 iteration_start_time_in_secs_ = start_time_in_secs_;
72 Init(options, parameters, solver_summary);
73 RETURN_IF_ERROR_AND_LOG(IterationZero());
74
75 // Create the TrustRegionStepEvaluator. The construction needs to be
76 // delayed to this point because we need the cost for the starting
77 // point to initialize the step evaluator.
78 step_evaluator_.reset(new TrustRegionStepEvaluator(
79 x_cost_,
80 options_.use_nonmonotonic_steps
81 ? options_.max_consecutive_nonmonotonic_steps
82 : 0));
83
84 while (FinalizeIterationAndCheckIfMinimizerCanContinue()) {
85 iteration_start_time_in_secs_ = WallTimeInSeconds();
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080086
87 const double previous_gradient_norm = iteration_summary_.gradient_norm;
88 const double previous_gradient_max_norm =
89 iteration_summary_.gradient_max_norm;
90
Austin Schuh70cc9552019-01-21 19:46:48 -080091 iteration_summary_ = IterationSummary();
92 iteration_summary_.iteration =
93 solver_summary->iterations.back().iteration + 1;
94
95 RETURN_IF_ERROR_AND_LOG(ComputeTrustRegionStep());
96 if (!iteration_summary_.step_is_valid) {
97 RETURN_IF_ERROR_AND_LOG(HandleInvalidStep());
98 continue;
99 }
100
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800101 if (options_.is_constrained &&
102 options_.max_num_line_search_step_size_iterations > 0) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800103 // Use a projected line search to enforce the bounds constraints
104 // and improve the quality of the step.
105 DoLineSearch(x_, gradient_, x_cost_, &delta_);
106 }
107
108 ComputeCandidatePointAndEvaluateCost();
109 DoInnerIterationsIfNeeded();
110
111 if (ParameterToleranceReached()) {
112 return;
113 }
114
115 if (FunctionToleranceReached()) {
116 return;
117 }
118
119 if (IsStepSuccessful()) {
120 RETURN_IF_ERROR_AND_LOG(HandleSuccessfulStep());
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800121 } else {
122 // Declare the step unsuccessful and inform the trust region strategy.
123 iteration_summary_.step_is_successful = false;
124 iteration_summary_.cost = candidate_cost_ + solver_summary_->fixed_cost;
Austin Schuh70cc9552019-01-21 19:46:48 -0800125
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800126 // When the step is unsuccessful, we do not compute the gradient
127 // (or update x), so we preserve its value from the last
128 // successful iteration.
129 iteration_summary_.gradient_norm = previous_gradient_norm;
130 iteration_summary_.gradient_max_norm = previous_gradient_max_norm;
131 strategy_->StepRejected(iteration_summary_.relative_decrease);
132 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800133 }
134}
135
136// Initialize the minimizer, allocate working space and set some of
137// the fields in the solver_summary.
138void TrustRegionMinimizer::Init(const Minimizer::Options& options,
139 double* parameters,
140 Solver::Summary* solver_summary) {
141 options_ = options;
142 sort(options_.trust_region_minimizer_iterations_to_dump.begin(),
143 options_.trust_region_minimizer_iterations_to_dump.end());
144
145 parameters_ = parameters;
146
147 solver_summary_ = solver_summary;
148 solver_summary_->termination_type = NO_CONVERGENCE;
149 solver_summary_->num_successful_steps = 0;
150 solver_summary_->num_unsuccessful_steps = 0;
151 solver_summary_->is_constrained = options.is_constrained;
152
153 CHECK(options_.evaluator != nullptr);
154 CHECK(options_.jacobian != nullptr);
155 CHECK(options_.trust_region_strategy != nullptr);
156 evaluator_ = options_.evaluator.get();
157 jacobian_ = options_.jacobian.get();
158 strategy_ = options_.trust_region_strategy.get();
159
160 is_not_silent_ = !options.is_silent;
161 inner_iterations_are_enabled_ =
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800162 options.inner_iteration_minimizer.get() != nullptr;
Austin Schuh70cc9552019-01-21 19:46:48 -0800163 inner_iterations_were_useful_ = false;
164
165 num_parameters_ = evaluator_->NumParameters();
166 num_effective_parameters_ = evaluator_->NumEffectiveParameters();
167 num_residuals_ = evaluator_->NumResiduals();
168 num_consecutive_invalid_steps_ = 0;
169
170 x_ = ConstVectorRef(parameters_, num_parameters_);
171 x_norm_ = x_.norm();
172 residuals_.resize(num_residuals_);
173 trust_region_step_.resize(num_effective_parameters_);
174 delta_.resize(num_effective_parameters_);
175 candidate_x_.resize(num_parameters_);
176 gradient_.resize(num_effective_parameters_);
177 model_residuals_.resize(num_residuals_);
178 negative_gradient_.resize(num_effective_parameters_);
179 projected_gradient_step_.resize(num_parameters_);
180
181 // By default scaling is one, if the user requests Jacobi scaling of
182 // the Jacobian, we will compute and overwrite this vector.
183 jacobian_scaling_ = Vector::Ones(num_effective_parameters_);
184
185 x_norm_ = -1; // Invalid value
186 x_cost_ = std::numeric_limits<double>::max();
187 minimum_cost_ = x_cost_;
188 model_cost_change_ = 0.0;
189}
190
191// 1. Project the initial solution onto the feasible set if needed.
192// 2. Compute the initial cost, jacobian & gradient.
193//
194// Return true if all computations can be performed successfully.
195bool TrustRegionMinimizer::IterationZero() {
196 iteration_summary_ = IterationSummary();
197 iteration_summary_.iteration = 0;
198 iteration_summary_.step_is_valid = false;
199 iteration_summary_.step_is_successful = false;
200 iteration_summary_.cost_change = 0.0;
201 iteration_summary_.gradient_max_norm = 0.0;
202 iteration_summary_.gradient_norm = 0.0;
203 iteration_summary_.step_norm = 0.0;
204 iteration_summary_.relative_decrease = 0.0;
205 iteration_summary_.eta = options_.eta;
206 iteration_summary_.linear_solver_iterations = 0;
207 iteration_summary_.step_solver_time_in_seconds = 0;
208
209 if (options_.is_constrained) {
210 delta_.setZero();
211 if (!evaluator_->Plus(x_.data(), delta_.data(), candidate_x_.data())) {
212 solver_summary_->message =
213 "Unable to project initial point onto the feasible set.";
214 solver_summary_->termination_type = FAILURE;
215 return false;
216 }
217
218 x_ = candidate_x_;
219 x_norm_ = x_.norm();
220 }
221
222 if (!EvaluateGradientAndJacobian(/*new_evaluation_point=*/true)) {
223 return false;
224 }
225
226 solver_summary_->initial_cost = x_cost_ + solver_summary_->fixed_cost;
227 iteration_summary_.step_is_valid = true;
228 iteration_summary_.step_is_successful = true;
229 return true;
230}
231
232// For the current x_, compute
233//
234// 1. Cost
235// 2. Jacobian
236// 3. Gradient
237// 4. Scale the Jacobian if needed (and compute the scaling if we are
238// in iteration zero).
239// 5. Compute the 2 and max norm of the gradient.
240//
241// Returns true if all computations could be performed
242// successfully. Any failures are considered fatal and the
243// Solver::Summary is updated to indicate this.
244bool TrustRegionMinimizer::EvaluateGradientAndJacobian(
245 bool new_evaluation_point) {
246 Evaluator::EvaluateOptions evaluate_options;
247 evaluate_options.new_evaluation_point = new_evaluation_point;
248 if (!evaluator_->Evaluate(evaluate_options,
249 x_.data(),
250 &x_cost_,
251 residuals_.data(),
252 gradient_.data(),
253 jacobian_)) {
254 solver_summary_->message = "Residual and Jacobian evaluation failed.";
255 solver_summary_->termination_type = FAILURE;
256 return false;
257 }
258
259 iteration_summary_.cost = x_cost_ + solver_summary_->fixed_cost;
260
261 if (options_.jacobi_scaling) {
262 if (iteration_summary_.iteration == 0) {
263 // Compute a scaling vector that is used to improve the
264 // conditioning of the Jacobian.
265 //
266 // jacobian_scaling_ = diag(J'J)^{-1}
267 jacobian_->SquaredColumnNorm(jacobian_scaling_.data());
268 for (int i = 0; i < jacobian_->num_cols(); ++i) {
269 // Add one to the denominator to prevent division by zero.
270 jacobian_scaling_[i] = 1.0 / (1.0 + sqrt(jacobian_scaling_[i]));
271 }
272 }
273
274 // jacobian = jacobian * diag(J'J) ^{-1}
275 jacobian_->ScaleColumns(jacobian_scaling_.data());
276 }
277
278 // The gradient exists in the local tangent space. To account for
279 // the bounds constraints correctly, instead of just computing the
280 // norm of the gradient vector, we compute
281 //
282 // |Plus(x, -gradient) - x|
283 //
284 // Where the Plus operator lifts the negative gradient to the
285 // ambient space, adds it to x and projects it on the hypercube
286 // defined by the bounds.
287 negative_gradient_ = -gradient_;
288 if (!evaluator_->Plus(x_.data(),
289 negative_gradient_.data(),
290 projected_gradient_step_.data())) {
291 solver_summary_->message =
292 "projected_gradient_step = Plus(x, -gradient) failed.";
293 solver_summary_->termination_type = FAILURE;
294 return false;
295 }
296
297 iteration_summary_.gradient_max_norm =
298 (x_ - projected_gradient_step_).lpNorm<Eigen::Infinity>();
299 iteration_summary_.gradient_norm = (x_ - projected_gradient_step_).norm();
300 return true;
301}
302
303// 1. Add the final timing information to the iteration summary.
304// 2. Run the callbacks
305// 3. Check for termination based on
306// a. Run time
307// b. Iteration count
308// c. Max norm of the gradient
309// d. Size of the trust region radius.
310//
311// Returns true if user did not terminate the solver and none of these
312// termination criterion are met.
313bool TrustRegionMinimizer::FinalizeIterationAndCheckIfMinimizerCanContinue() {
314 if (iteration_summary_.step_is_successful) {
315 ++solver_summary_->num_successful_steps;
316 if (x_cost_ < minimum_cost_) {
317 minimum_cost_ = x_cost_;
318 VectorRef(parameters_, num_parameters_) = x_;
319 iteration_summary_.step_is_nonmonotonic = false;
320 } else {
321 iteration_summary_.step_is_nonmonotonic = true;
322 }
323 } else {
324 ++solver_summary_->num_unsuccessful_steps;
325 }
326
327 iteration_summary_.trust_region_radius = strategy_->Radius();
328 iteration_summary_.iteration_time_in_seconds =
329 WallTimeInSeconds() - iteration_start_time_in_secs_;
330 iteration_summary_.cumulative_time_in_seconds =
331 WallTimeInSeconds() - start_time_in_secs_ +
332 solver_summary_->preprocessor_time_in_seconds;
333
334 solver_summary_->iterations.push_back(iteration_summary_);
335
336 if (!RunCallbacks(options_, iteration_summary_, solver_summary_)) {
337 return false;
338 }
339
340 if (MaxSolverTimeReached()) {
341 return false;
342 }
343
344 if (MaxSolverIterationsReached()) {
345 return false;
346 }
347
348 if (GradientToleranceReached()) {
349 return false;
350 }
351
352 if (MinTrustRegionRadiusReached()) {
353 return false;
354 }
355
356 return true;
357}
358
359// Compute the trust region step using the TrustRegionStrategy chosen
360// by the user.
361//
362// If the strategy returns with LINEAR_SOLVER_FATAL_ERROR, which
363// indicates an unrecoverable error, return false. This is the only
364// condition that returns false.
365//
366// If the strategy returns with LINEAR_SOLVER_FAILURE, which indicates
367// a numerical failure that could be recovered from by retrying
368// (e.g. by increasing the strength of the regularization), we set
369// iteration_summary_.step_is_valid to false and return true.
370//
371// In all other cases, we compute the decrease in the trust region
372// model problem. In exact arithmetic, this should always be
373// positive, but due to numerical problems in the TrustRegionStrategy
374// or round off error when computing the decrease it may be
375// negative. In which case again, we set
376// iteration_summary_.step_is_valid to false.
377bool TrustRegionMinimizer::ComputeTrustRegionStep() {
378 const double strategy_start_time = WallTimeInSeconds();
379 iteration_summary_.step_is_valid = false;
380 TrustRegionStrategy::PerSolveOptions per_solve_options;
381 per_solve_options.eta = options_.eta;
382 if (find(options_.trust_region_minimizer_iterations_to_dump.begin(),
383 options_.trust_region_minimizer_iterations_to_dump.end(),
384 iteration_summary_.iteration) !=
385 options_.trust_region_minimizer_iterations_to_dump.end()) {
386 per_solve_options.dump_format_type =
387 options_.trust_region_problem_dump_format_type;
388 per_solve_options.dump_filename_base =
389 JoinPath(options_.trust_region_problem_dump_directory,
390 StringPrintf("ceres_solver_iteration_%03d",
391 iteration_summary_.iteration));
392 }
393
394 TrustRegionStrategy::Summary strategy_summary =
395 strategy_->ComputeStep(per_solve_options,
396 jacobian_,
397 residuals_.data(),
398 trust_region_step_.data());
399
400 if (strategy_summary.termination_type == LINEAR_SOLVER_FATAL_ERROR) {
401 solver_summary_->message =
402 "Linear solver failed due to unrecoverable "
403 "non-numeric causes. Please see the error log for clues. ";
404 solver_summary_->termination_type = FAILURE;
405 return false;
406 }
407
408 iteration_summary_.step_solver_time_in_seconds =
409 WallTimeInSeconds() - strategy_start_time;
410 iteration_summary_.linear_solver_iterations = strategy_summary.num_iterations;
411
412 if (strategy_summary.termination_type == LINEAR_SOLVER_FAILURE) {
413 return true;
414 }
415
416 // new_model_cost
417 // = 1/2 [f + J * step]^2
418 // = 1/2 [ f'f + 2f'J * step + step' * J' * J * step ]
419 // model_cost_change
420 // = cost - new_model_cost
421 // = f'f/2 - 1/2 [ f'f + 2f'J * step + step' * J' * J * step]
422 // = -f'J * step - step' * J' * J * step / 2
423 // = -(J * step)'(f + J * step / 2)
424 model_residuals_.setZero();
425 jacobian_->RightMultiply(trust_region_step_.data(), model_residuals_.data());
426 model_cost_change_ =
427 -model_residuals_.dot(residuals_ + model_residuals_ / 2.0);
428
429 // TODO(sameeragarwal)
430 //
431 // 1. What happens if model_cost_change_ = 0
432 // 2. What happens if -epsilon <= model_cost_change_ < 0 for some
433 // small epsilon due to round off error.
434 iteration_summary_.step_is_valid = (model_cost_change_ > 0.0);
435 if (iteration_summary_.step_is_valid) {
436 // Undo the Jacobian column scaling.
437 delta_ = (trust_region_step_.array() * jacobian_scaling_.array()).matrix();
438 num_consecutive_invalid_steps_ = 0;
439 }
440
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800441 if (is_not_silent_ && !iteration_summary_.step_is_valid) {
442 VLOG(1) << "Invalid step: current_cost: " << x_cost_
443 << " absolute model cost change: " << model_cost_change_
444 << " relative model cost change: "
445 << (model_cost_change_ / x_cost_);
446 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800447 return true;
448}
449
450// Invalid steps can happen due to a number of reasons, and we allow a
451// limited number of consecutive failures, and return false if this
452// limit is exceeded.
453bool TrustRegionMinimizer::HandleInvalidStep() {
454 // TODO(sameeragarwal): Should we be returning FAILURE or
455 // NO_CONVERGENCE? The solution value is still usable in many cases,
456 // it is not clear if we should declare the solver a failure
457 // entirely. For example the case where model_cost_change ~ 0.0, but
458 // just slightly negative.
459 if (++num_consecutive_invalid_steps_ >=
460 options_.max_num_consecutive_invalid_steps) {
461 solver_summary_->message = StringPrintf(
462 "Number of consecutive invalid steps more "
463 "than Solver::Options::max_num_consecutive_invalid_steps: %d",
464 options_.max_num_consecutive_invalid_steps);
465 solver_summary_->termination_type = FAILURE;
466 return false;
467 }
468
469 strategy_->StepIsInvalid();
470
471 // We are going to try and reduce the trust region radius and
472 // solve again. To do this, we are going to treat this iteration
473 // as an unsuccessful iteration. Since the various callbacks are
474 // still executed, we are going to fill the iteration summary
475 // with data that assumes a step of length zero and no progress.
476 iteration_summary_.cost = x_cost_ + solver_summary_->fixed_cost;
477 iteration_summary_.cost_change = 0.0;
478 iteration_summary_.gradient_max_norm =
479 solver_summary_->iterations.back().gradient_max_norm;
480 iteration_summary_.gradient_norm =
481 solver_summary_->iterations.back().gradient_norm;
482 iteration_summary_.step_norm = 0.0;
483 iteration_summary_.relative_decrease = 0.0;
484 iteration_summary_.eta = options_.eta;
485 return true;
486}
487
488// Use the supplied coordinate descent minimizer to perform inner
489// iterations and compute the improvement due to it. Returns the cost
490// after performing the inner iterations.
491//
492// The optimization is performed with candidate_x_ as the starting
493// point, and if the optimization is successful, candidate_x_ will be
494// updated with the optimized parameters.
495void TrustRegionMinimizer::DoInnerIterationsIfNeeded() {
496 inner_iterations_were_useful_ = false;
497 if (!inner_iterations_are_enabled_ ||
498 candidate_cost_ >= std::numeric_limits<double>::max()) {
499 return;
500 }
501
502 double inner_iteration_start_time = WallTimeInSeconds();
503 ++solver_summary_->num_inner_iteration_steps;
504 inner_iteration_x_ = candidate_x_;
505 Solver::Summary inner_iteration_summary;
506 options_.inner_iteration_minimizer->Minimize(
507 options_, inner_iteration_x_.data(), &inner_iteration_summary);
508 double inner_iteration_cost;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800509 if (!evaluator_->Evaluate(inner_iteration_x_.data(),
510 &inner_iteration_cost,
511 nullptr,
512 nullptr,
513 nullptr)) {
514 if (is_not_silent_) {
515 VLOG(2) << "Inner iteration failed.";
516 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800517 return;
518 }
519
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800520 if (is_not_silent_) {
521 VLOG(2) << "Inner iteration succeeded; Current cost: " << x_cost_
522 << " Trust region step cost: " << candidate_cost_
523 << " Inner iteration cost: " << inner_iteration_cost;
524 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800525 candidate_x_ = inner_iteration_x_;
526
527 // Normally, the quality of a trust region step is measured by
528 // the ratio
529 //
530 // cost_change
531 // r = -----------------
532 // model_cost_change
533 //
534 // All the change in the nonlinear objective is due to the trust
535 // region step so this ratio is a good measure of the quality of
536 // the trust region radius. However, when inner iterations are
537 // being used, cost_change includes the contribution of the
538 // inner iterations and its not fair to credit it all to the
539 // trust region algorithm. So we change the ratio to be
540 //
541 // cost_change
542 // r = ------------------------------------------------
543 // (model_cost_change + inner_iteration_cost_change)
544 //
545 // Practically we do this by increasing model_cost_change by
546 // inner_iteration_cost_change.
547
548 const double inner_iteration_cost_change =
549 candidate_cost_ - inner_iteration_cost;
550 model_cost_change_ += inner_iteration_cost_change;
551 inner_iterations_were_useful_ = inner_iteration_cost < x_cost_;
552 const double inner_iteration_relative_progress =
553 1.0 - inner_iteration_cost / candidate_cost_;
554
555 // Disable inner iterations once the relative improvement
556 // drops below tolerance.
557 inner_iterations_are_enabled_ =
558 (inner_iteration_relative_progress > options_.inner_iteration_tolerance);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800559 if (is_not_silent_ && !inner_iterations_are_enabled_) {
560 VLOG(2) << "Disabling inner iterations. Progress : "
561 << inner_iteration_relative_progress;
562 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800563 candidate_cost_ = inner_iteration_cost;
564
565 solver_summary_->inner_iteration_time_in_seconds +=
566 WallTimeInSeconds() - inner_iteration_start_time;
567}
568
569// Perform a projected line search to improve the objective function
570// value along delta.
571//
572// TODO(sameeragarwal): The current implementation does not do
573// anything illegal but is incorrect and not terribly effective.
574//
575// https://github.com/ceres-solver/ceres-solver/issues/187
576void TrustRegionMinimizer::DoLineSearch(const Vector& x,
577 const Vector& gradient,
578 const double cost,
579 Vector* delta) {
580 LineSearchFunction line_search_function(evaluator_);
581
582 LineSearch::Options line_search_options;
583 line_search_options.is_silent = true;
584 line_search_options.interpolation_type =
585 options_.line_search_interpolation_type;
586 line_search_options.min_step_size = options_.min_line_search_step_size;
587 line_search_options.sufficient_decrease =
588 options_.line_search_sufficient_function_decrease;
589 line_search_options.max_step_contraction =
590 options_.max_line_search_step_contraction;
591 line_search_options.min_step_contraction =
592 options_.min_line_search_step_contraction;
593 line_search_options.max_num_iterations =
594 options_.max_num_line_search_step_size_iterations;
595 line_search_options.sufficient_curvature_decrease =
596 options_.line_search_sufficient_curvature_decrease;
597 line_search_options.max_step_expansion =
598 options_.max_line_search_step_expansion;
599 line_search_options.function = &line_search_function;
600
601 std::string message;
602 std::unique_ptr<LineSearch> line_search(
603 LineSearch::Create(ceres::ARMIJO, line_search_options, &message));
604 LineSearch::Summary line_search_summary;
605 line_search_function.Init(x, *delta);
606 line_search->Search(1.0, cost, gradient.dot(*delta), &line_search_summary);
607
608 solver_summary_->num_line_search_steps += line_search_summary.num_iterations;
609 solver_summary_->line_search_cost_evaluation_time_in_seconds +=
610 line_search_summary.cost_evaluation_time_in_seconds;
611 solver_summary_->line_search_gradient_evaluation_time_in_seconds +=
612 line_search_summary.gradient_evaluation_time_in_seconds;
613 solver_summary_->line_search_polynomial_minimization_time_in_seconds +=
614 line_search_summary.polynomial_minimization_time_in_seconds;
615 solver_summary_->line_search_total_time_in_seconds +=
616 line_search_summary.total_time_in_seconds;
617
618 if (line_search_summary.success) {
619 *delta *= line_search_summary.optimal_point.x;
620 }
621}
622
623// Check if the maximum amount of time allowed by the user for the
624// solver has been exceeded, and if so return false after updating
625// Solver::Summary::message.
626bool TrustRegionMinimizer::MaxSolverTimeReached() {
627 const double total_solver_time =
628 WallTimeInSeconds() - start_time_in_secs_ +
629 solver_summary_->preprocessor_time_in_seconds;
630 if (total_solver_time < options_.max_solver_time_in_seconds) {
631 return false;
632 }
633
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800634 solver_summary_->message = StringPrintf(
635 "Maximum solver time reached. "
636 "Total solver time: %e >= %e.",
637 total_solver_time,
638 options_.max_solver_time_in_seconds);
Austin Schuh70cc9552019-01-21 19:46:48 -0800639 solver_summary_->termination_type = NO_CONVERGENCE;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800640 if (is_not_silent_) {
641 VLOG(1) << "Terminating: " << solver_summary_->message;
642 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800643 return true;
644}
645
646// Check if the maximum number of iterations allowed by the user for
647// the solver has been exceeded, and if so return false after updating
648// Solver::Summary::message.
649bool TrustRegionMinimizer::MaxSolverIterationsReached() {
650 if (iteration_summary_.iteration < options_.max_num_iterations) {
651 return false;
652 }
653
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800654 solver_summary_->message = StringPrintf(
655 "Maximum number of iterations reached. "
656 "Number of iterations: %d.",
657 iteration_summary_.iteration);
Austin Schuh70cc9552019-01-21 19:46:48 -0800658
659 solver_summary_->termination_type = NO_CONVERGENCE;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800660 if (is_not_silent_) {
661 VLOG(1) << "Terminating: " << solver_summary_->message;
662 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800663 return true;
664}
665
666// Check convergence based on the max norm of the gradient (only for
667// iterations where the step was declared successful).
668bool TrustRegionMinimizer::GradientToleranceReached() {
669 if (!iteration_summary_.step_is_successful ||
670 iteration_summary_.gradient_max_norm > options_.gradient_tolerance) {
671 return false;
672 }
673
674 solver_summary_->message = StringPrintf(
675 "Gradient tolerance reached. "
676 "Gradient max norm: %e <= %e",
677 iteration_summary_.gradient_max_norm,
678 options_.gradient_tolerance);
679 solver_summary_->termination_type = CONVERGENCE;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800680 if (is_not_silent_) {
681 VLOG(1) << "Terminating: " << solver_summary_->message;
682 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800683 return true;
684}
685
686// Check convergence based the size of the trust region radius.
687bool TrustRegionMinimizer::MinTrustRegionRadiusReached() {
688 if (iteration_summary_.trust_region_radius >
689 options_.min_trust_region_radius) {
690 return false;
691 }
692
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800693 solver_summary_->message = StringPrintf(
694 "Minimum trust region radius reached. "
695 "Trust region radius: %e <= %e",
696 iteration_summary_.trust_region_radius,
697 options_.min_trust_region_radius);
Austin Schuh70cc9552019-01-21 19:46:48 -0800698 solver_summary_->termination_type = CONVERGENCE;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800699 if (is_not_silent_) {
700 VLOG(1) << "Terminating: " << solver_summary_->message;
701 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800702 return true;
703}
704
705// Solver::Options::parameter_tolerance based convergence check.
706bool TrustRegionMinimizer::ParameterToleranceReached() {
707 // Compute the norm of the step in the ambient space.
708 iteration_summary_.step_norm = (x_ - candidate_x_).norm();
709 const double step_size_tolerance =
710 options_.parameter_tolerance * (x_norm_ + options_.parameter_tolerance);
711
712 if (iteration_summary_.step_norm > step_size_tolerance) {
713 return false;
714 }
715
716 solver_summary_->message = StringPrintf(
717 "Parameter tolerance reached. "
718 "Relative step_norm: %e <= %e.",
719 (iteration_summary_.step_norm / (x_norm_ + options_.parameter_tolerance)),
720 options_.parameter_tolerance);
721 solver_summary_->termination_type = CONVERGENCE;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800722 if (is_not_silent_) {
723 VLOG(1) << "Terminating: " << solver_summary_->message;
724 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800725 return true;
726}
727
728// Solver::Options::function_tolerance based convergence check.
729bool TrustRegionMinimizer::FunctionToleranceReached() {
730 iteration_summary_.cost_change = x_cost_ - candidate_cost_;
731 const double absolute_function_tolerance =
732 options_.function_tolerance * x_cost_;
733
734 if (fabs(iteration_summary_.cost_change) > absolute_function_tolerance) {
735 return false;
736 }
737
738 solver_summary_->message = StringPrintf(
739 "Function tolerance reached. "
740 "|cost_change|/cost: %e <= %e",
741 fabs(iteration_summary_.cost_change) / x_cost_,
742 options_.function_tolerance);
743 solver_summary_->termination_type = CONVERGENCE;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800744 if (is_not_silent_) {
745 VLOG(1) << "Terminating: " << solver_summary_->message;
746 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800747 return true;
748}
749
750// Compute candidate_x_ = Plus(x_, delta_)
751// Evaluate the cost of candidate_x_ as candidate_cost_.
752//
753// Failure to compute the step or the cost mean that candidate_cost_
754// is set to std::numeric_limits<double>::max(). Unlike
755// EvaluateGradientAndJacobian, failure in this function is not fatal
756// as we are only computing and evaluating a candidate point, and if
757// for some reason we are unable to evaluate it, we consider it to be
758// a point with very high cost. This allows the user to deal with edge
759// cases/constraints as part of the LocalParameterization and
760// CostFunction objects.
761void TrustRegionMinimizer::ComputeCandidatePointAndEvaluateCost() {
762 if (!evaluator_->Plus(x_.data(), delta_.data(), candidate_x_.data())) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800763 if (is_not_silent_) {
764 LOG(WARNING) << "x_plus_delta = Plus(x, delta) failed. "
765 << "Treating it as a step with infinite cost";
766 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800767 candidate_cost_ = std::numeric_limits<double>::max();
768 return;
769 }
770
771 if (!evaluator_->Evaluate(
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800772 candidate_x_.data(), &candidate_cost_, nullptr, nullptr, nullptr)) {
773 if (is_not_silent_) {
774 LOG(WARNING) << "Step failed to evaluate. "
775 << "Treating it as a step with infinite cost";
776 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800777 candidate_cost_ = std::numeric_limits<double>::max();
778 }
779}
780
781bool TrustRegionMinimizer::IsStepSuccessful() {
782 iteration_summary_.relative_decrease =
783 step_evaluator_->StepQuality(candidate_cost_, model_cost_change_);
784
785 // In most cases, boosting the model_cost_change by the
786 // improvement caused by the inner iterations is fine, but it can
787 // be the case that the original trust region step was so bad that
788 // the resulting improvement in the cost was negative, and the
789 // change caused by the inner iterations was large enough to
790 // improve the step, but also to make relative decrease quite
791 // small.
792 //
793 // This can cause the trust region loop to reject this step. To
794 // get around this, we explicitly check if the inner iterations
795 // led to a net decrease in the objective function value. If
796 // they did, we accept the step even if the trust region ratio
797 // is small.
798 //
799 // Notice that we do not just check that cost_change is positive
800 // which is a weaker condition and would render the
801 // min_relative_decrease threshold useless. Instead, we keep
802 // track of inner_iterations_were_useful, which is true only
803 // when inner iterations lead to a net decrease in the cost.
804 return (inner_iterations_were_useful_ ||
805 iteration_summary_.relative_decrease >
806 options_.min_relative_decrease);
807}
808
809// Declare the step successful, move to candidate_x, update the
810// derivatives and let the trust region strategy and the step
811// evaluator know that the step has been accepted.
812bool TrustRegionMinimizer::HandleSuccessfulStep() {
813 x_ = candidate_x_;
814 x_norm_ = x_.norm();
815
816 // Since the step was successful, this point has already had the residual
817 // evaluated (but not the jacobian). So indicate that to the evaluator.
818 if (!EvaluateGradientAndJacobian(/*new_evaluation_point=*/false)) {
819 return false;
820 }
821
822 iteration_summary_.step_is_successful = true;
823 strategy_->StepAccepted(iteration_summary_.relative_decrease);
824 step_evaluator_->StepAccepted(candidate_cost_, model_cost_change_);
825 return true;
826}
827
Austin Schuh70cc9552019-01-21 19:46:48 -0800828} // namespace internal
829} // namespace ceres