dogleg_strategy.cc 27 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718
  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/dogleg_strategy.h"
  31. #include <algorithm>
  32. #include <cmath>
  33. #include "Eigen/Dense"
  34. #include "ceres/array_utils.h"
  35. #include "ceres/internal/eigen.h"
  36. #include "ceres/linear_least_squares_problems.h"
  37. #include "ceres/linear_solver.h"
  38. #include "ceres/polynomial.h"
  39. #include "ceres/sparse_matrix.h"
  40. #include "ceres/trust_region_strategy.h"
  41. #include "ceres/types.h"
  42. #include "glog/logging.h"
  43. namespace ceres::internal {
  44. namespace {
  45. const double kMaxMu = 1.0;
  46. const double kMinMu = 1e-8;
  47. } // namespace
  48. DoglegStrategy::DoglegStrategy(const TrustRegionStrategy::Options& options)
  49. : linear_solver_(options.linear_solver),
  50. radius_(options.initial_radius),
  51. max_radius_(options.max_radius),
  52. min_diagonal_(options.min_lm_diagonal),
  53. max_diagonal_(options.max_lm_diagonal),
  54. mu_(kMinMu),
  55. min_mu_(kMinMu),
  56. max_mu_(kMaxMu),
  57. mu_increase_factor_(10.0),
  58. increase_threshold_(0.75),
  59. decrease_threshold_(0.25),
  60. dogleg_step_norm_(0.0),
  61. reuse_(false),
  62. dogleg_type_(options.dogleg_type) {
  63. CHECK(linear_solver_ != nullptr);
  64. CHECK_GT(min_diagonal_, 0.0);
  65. CHECK_LE(min_diagonal_, max_diagonal_);
  66. CHECK_GT(max_radius_, 0.0);
  67. }
  68. // If the reuse_ flag is not set, then the Cauchy point (scaled
  69. // gradient) and the new Gauss-Newton step are computed from
  70. // scratch. The Dogleg step is then computed as interpolation of these
  71. // two vectors.
  72. TrustRegionStrategy::Summary DoglegStrategy::ComputeStep(
  73. const TrustRegionStrategy::PerSolveOptions& per_solve_options,
  74. SparseMatrix* jacobian,
  75. const double* residuals,
  76. double* step) {
  77. CHECK(jacobian != nullptr);
  78. CHECK(residuals != nullptr);
  79. CHECK(step != nullptr);
  80. const int n = jacobian->num_cols();
  81. if (reuse_) {
  82. // Gauss-Newton and gradient vectors are always available, only a
  83. // new interpolant need to be computed. For the subspace case,
  84. // the subspace and the two-dimensional model are also still valid.
  85. switch (dogleg_type_) {
  86. case TRADITIONAL_DOGLEG:
  87. ComputeTraditionalDoglegStep(step);
  88. break;
  89. case SUBSPACE_DOGLEG:
  90. ComputeSubspaceDoglegStep(step);
  91. break;
  92. }
  93. TrustRegionStrategy::Summary summary;
  94. summary.num_iterations = 0;
  95. summary.termination_type = LinearSolverTerminationType::SUCCESS;
  96. return summary;
  97. }
  98. reuse_ = true;
  99. // Check that we have the storage needed to hold the various
  100. // temporary vectors.
  101. if (diagonal_.rows() != n) {
  102. diagonal_.resize(n, 1);
  103. gradient_.resize(n, 1);
  104. gauss_newton_step_.resize(n, 1);
  105. }
  106. // Vector used to form the diagonal matrix that is used to
  107. // regularize the Gauss-Newton solve and that defines the
  108. // elliptical trust region
  109. //
  110. // || D * step || <= radius_ .
  111. //
  112. jacobian->SquaredColumnNorm(diagonal_.data());
  113. for (int i = 0; i < n; ++i) {
  114. diagonal_[i] =
  115. std::min(std::max(diagonal_[i], min_diagonal_), max_diagonal_);
  116. }
  117. diagonal_ = diagonal_.array().sqrt();
  118. ComputeGradient(jacobian, residuals);
  119. ComputeCauchyPoint(jacobian);
  120. LinearSolver::Summary linear_solver_summary =
  121. ComputeGaussNewtonStep(per_solve_options, jacobian, residuals);
  122. TrustRegionStrategy::Summary summary;
  123. summary.residual_norm = linear_solver_summary.residual_norm;
  124. summary.num_iterations = linear_solver_summary.num_iterations;
  125. summary.termination_type = linear_solver_summary.termination_type;
  126. if (linear_solver_summary.termination_type ==
  127. LinearSolverTerminationType::FATAL_ERROR) {
  128. return summary;
  129. }
  130. if (linear_solver_summary.termination_type !=
  131. LinearSolverTerminationType::FAILURE) {
  132. switch (dogleg_type_) {
  133. // Interpolate the Cauchy point and the Gauss-Newton step.
  134. case TRADITIONAL_DOGLEG:
  135. ComputeTraditionalDoglegStep(step);
  136. break;
  137. // Find the minimum in the subspace defined by the
  138. // Cauchy point and the (Gauss-)Newton step.
  139. case SUBSPACE_DOGLEG:
  140. if (!ComputeSubspaceModel(jacobian)) {
  141. summary.termination_type = LinearSolverTerminationType::FAILURE;
  142. break;
  143. }
  144. ComputeSubspaceDoglegStep(step);
  145. break;
  146. }
  147. }
  148. return summary;
  149. }
  150. // The trust region is assumed to be elliptical with the
  151. // diagonal scaling matrix D defined by sqrt(diagonal_).
  152. // It is implemented by substituting step' = D * step.
  153. // The trust region for step' is spherical.
  154. // The gradient, the Gauss-Newton step, the Cauchy point,
  155. // and all calculations involving the Jacobian have to
  156. // be adjusted accordingly.
  157. void DoglegStrategy::ComputeGradient(SparseMatrix* jacobian,
  158. const double* residuals) {
  159. gradient_.setZero();
  160. jacobian->LeftMultiplyAndAccumulate(residuals, gradient_.data());
  161. gradient_.array() /= diagonal_.array();
  162. }
  163. // The Cauchy point is the global minimizer of the quadratic model
  164. // along the one-dimensional subspace spanned by the gradient.
  165. void DoglegStrategy::ComputeCauchyPoint(SparseMatrix* jacobian) {
  166. // alpha * -gradient is the Cauchy point.
  167. Vector Jg(jacobian->num_rows());
  168. Jg.setZero();
  169. // The Jacobian is scaled implicitly by computing J * (D^-1 * (D^-1 * g))
  170. // instead of (J * D^-1) * (D^-1 * g).
  171. Vector scaled_gradient = (gradient_.array() / diagonal_.array()).matrix();
  172. jacobian->RightMultiplyAndAccumulate(scaled_gradient.data(), Jg.data());
  173. alpha_ = gradient_.squaredNorm() / Jg.squaredNorm();
  174. }
  175. // The dogleg step is defined as the intersection of the trust region
  176. // boundary with the piecewise linear path from the origin to the Cauchy
  177. // point and then from there to the Gauss-Newton point (global minimizer
  178. // of the model function). The Gauss-Newton point is taken if it lies
  179. // within the trust region.
  180. void DoglegStrategy::ComputeTraditionalDoglegStep(double* dogleg) {
  181. VectorRef dogleg_step(dogleg, gradient_.rows());
  182. // Case 1. The Gauss-Newton step lies inside the trust region, and
  183. // is therefore the optimal solution to the trust-region problem.
  184. const double gradient_norm = gradient_.norm();
  185. const double gauss_newton_norm = gauss_newton_step_.norm();
  186. if (gauss_newton_norm <= radius_) {
  187. dogleg_step = gauss_newton_step_;
  188. dogleg_step_norm_ = gauss_newton_norm;
  189. dogleg_step.array() /= diagonal_.array();
  190. VLOG(3) << "GaussNewton step size: " << dogleg_step_norm_
  191. << " radius: " << radius_;
  192. return;
  193. }
  194. // Case 2. The Cauchy point and the Gauss-Newton steps lie outside
  195. // the trust region. Rescale the Cauchy point to the trust region
  196. // and return.
  197. if (gradient_norm * alpha_ >= radius_) {
  198. dogleg_step = -(radius_ / gradient_norm) * gradient_;
  199. dogleg_step_norm_ = radius_;
  200. dogleg_step.array() /= diagonal_.array();
  201. VLOG(3) << "Cauchy step size: " << dogleg_step_norm_
  202. << " radius: " << radius_;
  203. return;
  204. }
  205. // Case 3. The Cauchy point is inside the trust region and the
  206. // Gauss-Newton step is outside. Compute the line joining the two
  207. // points and the point on it which intersects the trust region
  208. // boundary.
  209. // a = alpha * -gradient
  210. // b = gauss_newton_step
  211. const double b_dot_a = -alpha_ * gradient_.dot(gauss_newton_step_);
  212. const double a_squared_norm = pow(alpha_ * gradient_norm, 2.0);
  213. const double b_minus_a_squared_norm =
  214. a_squared_norm - 2 * b_dot_a + pow(gauss_newton_norm, 2);
  215. // c = a' (b - a)
  216. // = alpha * -gradient' gauss_newton_step - alpha^2 |gradient|^2
  217. const double c = b_dot_a - a_squared_norm;
  218. const double d = sqrt(c * c + b_minus_a_squared_norm *
  219. (pow(radius_, 2.0) - a_squared_norm));
  220. double beta = (c <= 0) ? (d - c) / b_minus_a_squared_norm
  221. : (radius_ * radius_ - a_squared_norm) / (d + c);
  222. dogleg_step =
  223. (-alpha_ * (1.0 - beta)) * gradient_ + beta * gauss_newton_step_;
  224. dogleg_step_norm_ = dogleg_step.norm();
  225. dogleg_step.array() /= diagonal_.array();
  226. VLOG(3) << "Dogleg step size: " << dogleg_step_norm_
  227. << " radius: " << radius_;
  228. }
  229. // The subspace method finds the minimum of the two-dimensional problem
  230. //
  231. // min. 1/2 x' B' H B x + g' B x
  232. // s.t. || B x ||^2 <= r^2
  233. //
  234. // where r is the trust region radius and B is the matrix with unit columns
  235. // spanning the subspace defined by the steepest descent and Newton direction.
  236. // This subspace by definition includes the Gauss-Newton point, which is
  237. // therefore taken if it lies within the trust region.
  238. void DoglegStrategy::ComputeSubspaceDoglegStep(double* dogleg) {
  239. VectorRef dogleg_step(dogleg, gradient_.rows());
  240. // The Gauss-Newton point is inside the trust region if |GN| <= radius_.
  241. // This test is valid even though radius_ is a length in the two-dimensional
  242. // subspace while gauss_newton_step_ is expressed in the (scaled)
  243. // higher dimensional original space. This is because
  244. //
  245. // 1. gauss_newton_step_ by definition lies in the subspace, and
  246. // 2. the subspace basis is orthonormal.
  247. //
  248. // As a consequence, the norm of the gauss_newton_step_ in the subspace is
  249. // the same as its norm in the original space.
  250. const double gauss_newton_norm = gauss_newton_step_.norm();
  251. if (gauss_newton_norm <= radius_) {
  252. dogleg_step = gauss_newton_step_;
  253. dogleg_step_norm_ = gauss_newton_norm;
  254. dogleg_step.array() /= diagonal_.array();
  255. VLOG(3) << "GaussNewton step size: " << dogleg_step_norm_
  256. << " radius: " << radius_;
  257. return;
  258. }
  259. // The optimum lies on the boundary of the trust region. The above problem
  260. // therefore becomes
  261. //
  262. // min. 1/2 x^T B^T H B x + g^T B x
  263. // s.t. || B x ||^2 = r^2
  264. //
  265. // Notice the equality in the constraint.
  266. //
  267. // This can be solved by forming the Lagrangian, solving for x(y), where
  268. // y is the Lagrange multiplier, using the gradient of the objective, and
  269. // putting x(y) back into the constraint. This results in a fourth order
  270. // polynomial in y, which can be solved using e.g. the companion matrix.
  271. // See the description of MakePolynomialForBoundaryConstrainedProblem for
  272. // details. The result is up to four real roots y*, not all of which
  273. // correspond to feasible points. The feasible points x(y*) have to be
  274. // tested for optimality.
  275. if (subspace_is_one_dimensional_) {
  276. // The subspace is one-dimensional, so both the gradient and
  277. // the Gauss-Newton step point towards the same direction.
  278. // In this case, we move along the gradient until we reach the trust
  279. // region boundary.
  280. dogleg_step = -(radius_ / gradient_.norm()) * gradient_;
  281. dogleg_step_norm_ = radius_;
  282. dogleg_step.array() /= diagonal_.array();
  283. VLOG(3) << "Dogleg subspace step size (1D): " << dogleg_step_norm_
  284. << " radius: " << radius_;
  285. return;
  286. }
  287. Vector2d minimum(0.0, 0.0);
  288. if (!FindMinimumOnTrustRegionBoundary(&minimum)) {
  289. // For the positive semi-definite case, a traditional dogleg step
  290. // is taken in this case.
  291. LOG(WARNING) << "Failed to compute polynomial roots. "
  292. << "Taking traditional dogleg step instead.";
  293. ComputeTraditionalDoglegStep(dogleg);
  294. return;
  295. }
  296. // Test first order optimality at the minimum.
  297. // The first order KKT conditions state that the minimum x*
  298. // has to satisfy either || x* ||^2 < r^2 (i.e. has to lie within
  299. // the trust region), or
  300. //
  301. // (B x* + g) + y x* = 0
  302. //
  303. // for some positive scalar y.
  304. // Here, as it is already known that the minimum lies on the boundary, the
  305. // latter condition is tested. To allow for small imprecisions, we test if
  306. // the angle between (B x* + g) and -x* is smaller than acos(0.99).
  307. // The exact value of the cosine is arbitrary but should be close to 1.
  308. //
  309. // This condition should not be violated. If it is, the minimum was not
  310. // correctly determined.
  311. const double kCosineThreshold = 0.99;
  312. const Vector2d grad_minimum = subspace_B_ * minimum + subspace_g_;
  313. const double cosine_angle =
  314. -minimum.dot(grad_minimum) / (minimum.norm() * grad_minimum.norm());
  315. if (cosine_angle < kCosineThreshold) {
  316. LOG(WARNING) << "First order optimality seems to be violated "
  317. << "in the subspace method!\n"
  318. << "Cosine of angle between x and B x + g is " << cosine_angle
  319. << ".\n"
  320. << "Taking a regular dogleg step instead.\n"
  321. << "Please consider filing a bug report if this "
  322. << "happens frequently or consistently.\n";
  323. ComputeTraditionalDoglegStep(dogleg);
  324. return;
  325. }
  326. // Create the full step from the optimal 2d solution.
  327. dogleg_step = subspace_basis_ * minimum;
  328. dogleg_step_norm_ = radius_;
  329. dogleg_step.array() /= diagonal_.array();
  330. VLOG(3) << "Dogleg subspace step size: " << dogleg_step_norm_
  331. << " radius: " << radius_;
  332. }
  333. // Build the polynomial that defines the optimal Lagrange multipliers.
  334. // Let the Lagrangian be
  335. //
  336. // L(x, y) = 0.5 x^T B x + x^T g + y (0.5 x^T x - 0.5 r^2). (1)
  337. //
  338. // Stationary points of the Lagrangian are given by
  339. //
  340. // 0 = d L(x, y) / dx = Bx + g + y x (2)
  341. // 0 = d L(x, y) / dy = 0.5 x^T x - 0.5 r^2 (3)
  342. //
  343. // For any given y, we can solve (2) for x as
  344. //
  345. // x(y) = -(B + y I)^-1 g . (4)
  346. //
  347. // As B + y I is 2x2, we form the inverse explicitly:
  348. //
  349. // (B + y I)^-1 = (1 / det(B + y I)) adj(B + y I) (5)
  350. //
  351. // where adj() denotes adjugation. This should be safe, as B is positive
  352. // semi-definite and y is necessarily positive, so (B + y I) is indeed
  353. // invertible.
  354. // Plugging (5) into (4) and the result into (3), then dividing by 0.5 we
  355. // obtain
  356. //
  357. // 0 = (1 / det(B + y I))^2 g^T adj(B + y I)^T adj(B + y I) g - r^2
  358. // (6)
  359. //
  360. // or
  361. //
  362. // det(B + y I)^2 r^2 = g^T adj(B + y I)^T adj(B + y I) g (7a)
  363. // = g^T adj(B)^T adj(B) g
  364. // + 2 y g^T adj(B)^T g + y^2 g^T g (7b)
  365. //
  366. // as
  367. //
  368. // adj(B + y I) = adj(B) + y I = adj(B)^T + y I . (8)
  369. //
  370. // The left hand side can be expressed explicitly using
  371. //
  372. // det(B + y I) = det(B) + y tr(B) + y^2 . (9)
  373. //
  374. // So (7) is a polynomial in y of degree four.
  375. // Bringing everything back to the left hand side, the coefficients can
  376. // be read off as
  377. //
  378. // y^4 r^2
  379. // + y^3 2 r^2 tr(B)
  380. // + y^2 (r^2 tr(B)^2 + 2 r^2 det(B) - g^T g)
  381. // + y^1 (2 r^2 det(B) tr(B) - 2 g^T adj(B)^T g)
  382. // + y^0 (r^2 det(B)^2 - g^T adj(B)^T adj(B) g)
  383. //
  384. Vector DoglegStrategy::MakePolynomialForBoundaryConstrainedProblem() const {
  385. const double detB = subspace_B_.determinant();
  386. const double trB = subspace_B_.trace();
  387. const double r2 = radius_ * radius_;
  388. Matrix2d B_adj;
  389. // clang-format off
  390. B_adj << subspace_B_(1, 1) , -subspace_B_(0, 1),
  391. -subspace_B_(1, 0) , subspace_B_(0, 0);
  392. // clang-format on
  393. Vector polynomial(5);
  394. polynomial(0) = r2;
  395. polynomial(1) = 2.0 * r2 * trB;
  396. polynomial(2) = r2 * (trB * trB + 2.0 * detB) - subspace_g_.squaredNorm();
  397. polynomial(3) =
  398. -2.0 * (subspace_g_.transpose() * B_adj * subspace_g_ - r2 * detB * trB);
  399. polynomial(4) = r2 * detB * detB - (B_adj * subspace_g_).squaredNorm();
  400. return polynomial;
  401. }
  402. // Given a Lagrange multiplier y that corresponds to a stationary point
  403. // of the Lagrangian L(x, y), compute the corresponding x from the
  404. // equation
  405. //
  406. // 0 = d L(x, y) / dx
  407. // = B * x + g + y * x
  408. // = (B + y * I) * x + g
  409. //
  410. DoglegStrategy::Vector2d DoglegStrategy::ComputeSubspaceStepFromRoot(
  411. double y) const {
  412. const Matrix2d B_i = subspace_B_ + y * Matrix2d::Identity();
  413. return -B_i.partialPivLu().solve(subspace_g_);
  414. }
  415. // This function evaluates the quadratic model at a point x in the
  416. // subspace spanned by subspace_basis_.
  417. double DoglegStrategy::EvaluateSubspaceModel(const Vector2d& x) const {
  418. return 0.5 * x.dot(subspace_B_ * x) + subspace_g_.dot(x);
  419. }
  420. // This function attempts to solve the boundary-constrained subspace problem
  421. //
  422. // min. 1/2 x^T B^T H B x + g^T B x
  423. // s.t. || B x ||^2 = r^2
  424. //
  425. // where B is an orthonormal subspace basis and r is the trust-region radius.
  426. //
  427. // This is done by finding the roots of a fourth degree polynomial. If the
  428. // root finding fails, the function returns false and minimum will be set
  429. // to (0, 0). If it succeeds, true is returned.
  430. //
  431. // In the failure case, another step should be taken, such as the traditional
  432. // dogleg step.
  433. bool DoglegStrategy::FindMinimumOnTrustRegionBoundary(Vector2d* minimum) const {
  434. CHECK(minimum != nullptr);
  435. // Return (0, 0) in all error cases.
  436. minimum->setZero();
  437. // Create the fourth-degree polynomial that is a necessary condition for
  438. // optimality.
  439. const Vector polynomial = MakePolynomialForBoundaryConstrainedProblem();
  440. // Find the real parts y_i of its roots (not only the real roots).
  441. Vector roots_real;
  442. if (!FindPolynomialRoots(polynomial, &roots_real, nullptr)) {
  443. // Failed to find the roots of the polynomial, i.e. the candidate
  444. // solutions of the constrained problem. Report this back to the caller.
  445. return false;
  446. }
  447. // For each root y, compute B x(y) and check for feasibility.
  448. // Notice that there should always be four roots, as the leading term of
  449. // the polynomial is r^2 and therefore non-zero. However, as some roots
  450. // may be complex, the real parts are not necessarily unique.
  451. double minimum_value = std::numeric_limits<double>::max();
  452. bool valid_root_found = false;
  453. for (int i = 0; i < roots_real.size(); ++i) {
  454. const Vector2d x_i = ComputeSubspaceStepFromRoot(roots_real(i));
  455. // Not all roots correspond to points on the trust region boundary.
  456. // There are at most four candidate solutions. As we are interested
  457. // in the minimum, it is safe to consider all of them after projecting
  458. // them onto the trust region boundary.
  459. if (x_i.norm() > 0) {
  460. const double f_i = EvaluateSubspaceModel((radius_ / x_i.norm()) * x_i);
  461. valid_root_found = true;
  462. if (f_i < minimum_value) {
  463. minimum_value = f_i;
  464. *minimum = x_i;
  465. }
  466. }
  467. }
  468. return valid_root_found;
  469. }
  470. LinearSolver::Summary DoglegStrategy::ComputeGaussNewtonStep(
  471. const PerSolveOptions& per_solve_options,
  472. SparseMatrix* jacobian,
  473. const double* residuals) {
  474. const int n = jacobian->num_cols();
  475. LinearSolver::Summary linear_solver_summary;
  476. linear_solver_summary.termination_type = LinearSolverTerminationType::FAILURE;
  477. // The Jacobian matrix is often quite poorly conditioned. Thus it is
  478. // necessary to add a diagonal matrix at the bottom to prevent the
  479. // linear solver from failing.
  480. //
  481. // We do this by computing the same diagonal matrix as the one used
  482. // by Levenberg-Marquardt (other choices are possible), and scaling
  483. // it by a small constant (independent of the trust region radius).
  484. //
  485. // If the solve fails, the multiplier to the diagonal is increased
  486. // up to max_mu_ by a factor of mu_increase_factor_ every time. If
  487. // the linear solver is still not successful, the strategy returns
  488. // with LinearSolverTerminationType::FAILURE.
  489. //
  490. // Next time when a new Gauss-Newton step is requested, the
  491. // multiplier starts out from the last successful solve.
  492. //
  493. // When a step is declared successful, the multiplier is decreased
  494. // by half of mu_increase_factor_.
  495. while (mu_ < max_mu_) {
  496. // Dogleg, as far as I (sameeragarwal) understand it, requires a
  497. // reasonably good estimate of the Gauss-Newton step. This means
  498. // that we need to solve the normal equations more or less
  499. // exactly. This is reflected in the values of the tolerances set
  500. // below.
  501. //
  502. // For now, this strategy should only be used with exact
  503. // factorization based solvers, for which these tolerances are
  504. // automatically satisfied.
  505. //
  506. // The right way to combine inexact solves with trust region
  507. // methods is to use Stiehaug's method.
  508. LinearSolver::PerSolveOptions solve_options;
  509. solve_options.q_tolerance = 0.0;
  510. solve_options.r_tolerance = 0.0;
  511. lm_diagonal_ = diagonal_ * std::sqrt(mu_);
  512. solve_options.D = lm_diagonal_.data();
  513. // As in the LevenbergMarquardtStrategy, solve Jy = r instead
  514. // of Jx = -r and later set x = -y to avoid having to modify
  515. // either jacobian or residuals.
  516. InvalidateArray(n, gauss_newton_step_.data());
  517. linear_solver_summary = linear_solver_->Solve(
  518. jacobian, residuals, solve_options, gauss_newton_step_.data());
  519. if (per_solve_options.dump_format_type == CONSOLE ||
  520. (per_solve_options.dump_format_type != CONSOLE &&
  521. !per_solve_options.dump_filename_base.empty())) {
  522. if (!DumpLinearLeastSquaresProblem(per_solve_options.dump_filename_base,
  523. per_solve_options.dump_format_type,
  524. jacobian,
  525. solve_options.D,
  526. residuals,
  527. gauss_newton_step_.data(),
  528. 0)) {
  529. LOG(ERROR) << "Unable to dump trust region problem."
  530. << " Filename base: "
  531. << per_solve_options.dump_filename_base;
  532. }
  533. }
  534. if (linear_solver_summary.termination_type ==
  535. LinearSolverTerminationType::FATAL_ERROR) {
  536. return linear_solver_summary;
  537. }
  538. if (linear_solver_summary.termination_type ==
  539. LinearSolverTerminationType::FAILURE ||
  540. !IsArrayValid(n, gauss_newton_step_.data())) {
  541. mu_ *= mu_increase_factor_;
  542. VLOG(2) << "Increasing mu " << mu_;
  543. linear_solver_summary.termination_type =
  544. LinearSolverTerminationType::FAILURE;
  545. continue;
  546. }
  547. break;
  548. }
  549. if (linear_solver_summary.termination_type !=
  550. LinearSolverTerminationType::FAILURE) {
  551. // The scaled Gauss-Newton step is D * GN:
  552. //
  553. // - (D^-1 J^T J D^-1)^-1 (D^-1 g)
  554. // = - D (J^T J)^-1 D D^-1 g
  555. // = D -(J^T J)^-1 g
  556. //
  557. gauss_newton_step_.array() *= -diagonal_.array();
  558. }
  559. return linear_solver_summary;
  560. }
  561. void DoglegStrategy::StepAccepted(double step_quality) {
  562. CHECK_GT(step_quality, 0.0);
  563. if (step_quality < decrease_threshold_) {
  564. radius_ *= 0.5;
  565. }
  566. if (step_quality > increase_threshold_) {
  567. radius_ = std::max(radius_, 3.0 * dogleg_step_norm_);
  568. }
  569. // Reduce the regularization multiplier, in the hope that whatever
  570. // was causing the rank deficiency has gone away and we can return
  571. // to doing a pure Gauss-Newton solve.
  572. mu_ = std::max(min_mu_, 2.0 * mu_ / mu_increase_factor_);
  573. reuse_ = false;
  574. }
  575. void DoglegStrategy::StepRejected(double /*step_quality*/) {
  576. radius_ *= 0.5;
  577. reuse_ = true;
  578. }
  579. void DoglegStrategy::StepIsInvalid() {
  580. mu_ *= mu_increase_factor_;
  581. reuse_ = false;
  582. }
  583. double DoglegStrategy::Radius() const { return radius_; }
  584. bool DoglegStrategy::ComputeSubspaceModel(SparseMatrix* jacobian) {
  585. // Compute an orthogonal basis for the subspace using QR decomposition.
  586. Matrix basis_vectors(jacobian->num_cols(), 2);
  587. basis_vectors.col(0) = gradient_;
  588. basis_vectors.col(1) = gauss_newton_step_;
  589. Eigen::ColPivHouseholderQR<Matrix> basis_qr(basis_vectors);
  590. switch (basis_qr.rank()) {
  591. case 0:
  592. // This should never happen, as it implies that both the gradient
  593. // and the Gauss-Newton step are zero. In this case, the minimizer should
  594. // have stopped due to the gradient being too small.
  595. LOG(ERROR) << "Rank of subspace basis is 0. "
  596. << "This means that the gradient at the current iterate is "
  597. << "zero but the optimization has not been terminated. "
  598. << "You may have found a bug in Ceres.";
  599. return false;
  600. case 1:
  601. // Gradient and Gauss-Newton step coincide, so we lie on one of the
  602. // major axes of the quadratic problem. In this case, we simply move
  603. // along the gradient until we reach the trust region boundary.
  604. subspace_is_one_dimensional_ = true;
  605. return true;
  606. case 2:
  607. subspace_is_one_dimensional_ = false;
  608. break;
  609. default:
  610. LOG(ERROR) << "Rank of the subspace basis matrix is reported to be "
  611. << "greater than 2. As the matrix contains only two "
  612. << "columns this cannot be true and is indicative of "
  613. << "a bug.";
  614. return false;
  615. }
  616. // The subspace is two-dimensional, so compute the subspace model.
  617. // Given the basis U, this is
  618. //
  619. // subspace_g_ = g_scaled^T U
  620. //
  621. // and
  622. //
  623. // subspace_B_ = U^T (J_scaled^T J_scaled) U
  624. //
  625. // As J_scaled = J * D^-1, the latter becomes
  626. //
  627. // subspace_B_ = ((U^T D^-1) J^T) (J (D^-1 U))
  628. // = (J (D^-1 U))^T (J (D^-1 U))
  629. subspace_basis_ =
  630. basis_qr.householderQ() * Matrix::Identity(jacobian->num_cols(), 2);
  631. subspace_g_ = subspace_basis_.transpose() * gradient_;
  632. Eigen::Matrix<double, 2, Eigen::Dynamic, Eigen::RowMajor> Jb(
  633. 2, jacobian->num_rows());
  634. Jb.setZero();
  635. Vector tmp;
  636. tmp = (subspace_basis_.col(0).array() / diagonal_.array()).matrix();
  637. jacobian->RightMultiplyAndAccumulate(tmp.data(), Jb.row(0).data());
  638. tmp = (subspace_basis_.col(1).array() / diagonal_.array()).matrix();
  639. jacobian->RightMultiplyAndAccumulate(tmp.data(), Jb.row(1).data());
  640. subspace_B_ = Jb * Jb.transpose();
  641. return true;
  642. }
  643. } // namespace ceres::internal