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