trust_region_minimizer.cc 32 KB

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