line_search_direction.cc 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369
  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/line_search_direction.h"
  31. #include <memory>
  32. #include "ceres/internal/eigen.h"
  33. #include "ceres/internal/export.h"
  34. #include "ceres/line_search_minimizer.h"
  35. #include "ceres/low_rank_inverse_hessian.h"
  36. #include "glog/logging.h"
  37. namespace ceres::internal {
  38. class CERES_NO_EXPORT SteepestDescent final : public LineSearchDirection {
  39. public:
  40. bool NextDirection(const LineSearchMinimizer::State& /*previous*/,
  41. const LineSearchMinimizer::State& current,
  42. Vector* search_direction) override {
  43. *search_direction = -current.gradient;
  44. return true;
  45. }
  46. };
  47. class CERES_NO_EXPORT NonlinearConjugateGradient final
  48. : public LineSearchDirection {
  49. public:
  50. NonlinearConjugateGradient(const NonlinearConjugateGradientType type,
  51. const double function_tolerance)
  52. : type_(type), function_tolerance_(function_tolerance) {}
  53. bool NextDirection(const LineSearchMinimizer::State& previous,
  54. const LineSearchMinimizer::State& current,
  55. Vector* search_direction) override {
  56. double beta = 0.0;
  57. Vector gradient_change;
  58. switch (type_) {
  59. case FLETCHER_REEVES:
  60. beta = current.gradient_squared_norm / previous.gradient_squared_norm;
  61. break;
  62. case POLAK_RIBIERE:
  63. gradient_change = current.gradient - previous.gradient;
  64. beta = (current.gradient.dot(gradient_change) /
  65. previous.gradient_squared_norm);
  66. break;
  67. case HESTENES_STIEFEL:
  68. gradient_change = current.gradient - previous.gradient;
  69. beta = (current.gradient.dot(gradient_change) /
  70. previous.search_direction.dot(gradient_change));
  71. break;
  72. default:
  73. LOG(FATAL) << "Unknown nonlinear conjugate gradient type: " << type_;
  74. }
  75. *search_direction = -current.gradient + beta * previous.search_direction;
  76. const double directional_derivative =
  77. current.gradient.dot(*search_direction);
  78. if (directional_derivative > -function_tolerance_) {
  79. LOG(WARNING) << "Restarting non-linear conjugate gradients: "
  80. << directional_derivative;
  81. *search_direction = -current.gradient;
  82. }
  83. return true;
  84. }
  85. private:
  86. const NonlinearConjugateGradientType type_;
  87. const double function_tolerance_;
  88. };
  89. class CERES_NO_EXPORT LBFGS final : public LineSearchDirection {
  90. public:
  91. LBFGS(const int num_parameters,
  92. const int max_lbfgs_rank,
  93. const bool use_approximate_eigenvalue_bfgs_scaling)
  94. : low_rank_inverse_hessian_(num_parameters,
  95. max_lbfgs_rank,
  96. use_approximate_eigenvalue_bfgs_scaling),
  97. is_positive_definite_(true) {}
  98. bool NextDirection(const LineSearchMinimizer::State& previous,
  99. const LineSearchMinimizer::State& current,
  100. Vector* search_direction) override {
  101. CHECK(is_positive_definite_)
  102. << "Ceres bug: NextDirection() called on L-BFGS after inverse Hessian "
  103. << "approximation has become indefinite, please contact the "
  104. << "developers!";
  105. low_rank_inverse_hessian_.Update(
  106. previous.search_direction * previous.step_size,
  107. current.gradient - previous.gradient);
  108. search_direction->setZero();
  109. low_rank_inverse_hessian_.RightMultiplyAndAccumulate(
  110. current.gradient.data(), search_direction->data());
  111. *search_direction *= -1.0;
  112. if (search_direction->dot(current.gradient) >= 0.0) {
  113. LOG(WARNING) << "Numerical failure in L-BFGS update: inverse Hessian "
  114. << "approximation is not positive definite, and thus "
  115. << "initial gradient for search direction is positive: "
  116. << search_direction->dot(current.gradient);
  117. is_positive_definite_ = false;
  118. return false;
  119. }
  120. return true;
  121. }
  122. private:
  123. LowRankInverseHessian low_rank_inverse_hessian_;
  124. bool is_positive_definite_;
  125. };
  126. class CERES_NO_EXPORT BFGS final : public LineSearchDirection {
  127. public:
  128. BFGS(const int num_parameters, const bool use_approximate_eigenvalue_scaling)
  129. : num_parameters_(num_parameters),
  130. use_approximate_eigenvalue_scaling_(use_approximate_eigenvalue_scaling),
  131. initialized_(false),
  132. is_positive_definite_(true) {
  133. if (num_parameters_ >= 1000) {
  134. LOG(WARNING) << "BFGS line search being created with: " << num_parameters_
  135. << " parameters, this will allocate a dense approximate "
  136. << "inverse Hessian of size: " << num_parameters_ << " x "
  137. << num_parameters_
  138. << ", consider using the L-BFGS memory-efficient line "
  139. << "search direction instead.";
  140. }
  141. // Construct inverse_hessian_ after logging warning about size s.t. if the
  142. // allocation crashes us, the log will highlight what the issue likely was.
  143. inverse_hessian_ = Matrix::Identity(num_parameters, num_parameters);
  144. }
  145. bool NextDirection(const LineSearchMinimizer::State& previous,
  146. const LineSearchMinimizer::State& current,
  147. Vector* search_direction) override {
  148. CHECK(is_positive_definite_)
  149. << "Ceres bug: NextDirection() called on BFGS after inverse Hessian "
  150. << "approximation has become indefinite, please contact the "
  151. << "developers!";
  152. const Vector delta_x = previous.search_direction * previous.step_size;
  153. const Vector delta_gradient = current.gradient - previous.gradient;
  154. const double delta_x_dot_delta_gradient = delta_x.dot(delta_gradient);
  155. // The (L)BFGS algorithm explicitly requires that the secant equation:
  156. //
  157. // B_{k+1} * s_k = y_k
  158. //
  159. // Is satisfied at each iteration, where B_{k+1} is the approximated
  160. // Hessian at the k+1-th iteration, s_k = (x_{k+1} - x_{k}) and
  161. // y_k = (grad_{k+1} - grad_{k}). As the approximated Hessian must be
  162. // positive definite, this is equivalent to the condition:
  163. //
  164. // s_k^T * y_k > 0 [s_k^T * B_{k+1} * s_k = s_k^T * y_k > 0]
  165. //
  166. // This condition would always be satisfied if the function was strictly
  167. // convex, alternatively, it is always satisfied provided that a Wolfe line
  168. // search is used (even if the function is not strictly convex). See [1]
  169. // (p138) for a proof.
  170. //
  171. // Although Ceres will always use a Wolfe line search when using (L)BFGS,
  172. // practical implementation considerations mean that the line search
  173. // may return a point that satisfies only the Armijo condition, and thus
  174. // could violate the Secant equation. As such, we will only use a step
  175. // to update the Hessian approximation if:
  176. //
  177. // s_k^T * y_k > tolerance
  178. //
  179. // It is important that tolerance is very small (and >=0), as otherwise we
  180. // might skip the update too often and fail to capture important curvature
  181. // information in the Hessian. For example going from 1e-10 -> 1e-14
  182. // improves the NIST benchmark score from 43/54 to 53/54.
  183. //
  184. // [1] Nocedal J, Wright S, Numerical Optimization, 2nd Ed. Springer, 1999.
  185. //
  186. // TODO(alexs.mac): Consider using Damped BFGS update instead of
  187. // skipping update.
  188. const double kBFGSSecantConditionHessianUpdateTolerance = 1e-14;
  189. if (delta_x_dot_delta_gradient <=
  190. kBFGSSecantConditionHessianUpdateTolerance) {
  191. VLOG(2) << "Skipping BFGS Update, delta_x_dot_delta_gradient too "
  192. << "small: " << delta_x_dot_delta_gradient
  193. << ", tolerance: " << kBFGSSecantConditionHessianUpdateTolerance
  194. << " (Secant condition).";
  195. } else {
  196. // Update dense inverse Hessian approximation.
  197. if (!initialized_ && use_approximate_eigenvalue_scaling_) {
  198. // Rescale the initial inverse Hessian approximation (H_0) to be
  199. // iteratively updated so that it is of similar 'size' to the true
  200. // inverse Hessian at the start point. As shown in [1]:
  201. //
  202. // \gamma = (delta_gradient_{0}' * delta_x_{0}) /
  203. // (delta_gradient_{0}' * delta_gradient_{0})
  204. //
  205. // Satisfies:
  206. //
  207. // (1 / \lambda_m) <= \gamma <= (1 / \lambda_1)
  208. //
  209. // Where \lambda_1 & \lambda_m are the smallest and largest eigenvalues
  210. // of the true initial Hessian (not the inverse) respectively. Thus,
  211. // \gamma is an approximate eigenvalue of the true inverse Hessian, and
  212. // choosing: H_0 = I * \gamma will yield a starting point that has a
  213. // similar scale to the true inverse Hessian. This technique is widely
  214. // reported to often improve convergence, however this is not
  215. // universally true, particularly if there are errors in the initial
  216. // gradients, or if there are significant differences in the sensitivity
  217. // of the problem to the parameters (i.e. the range of the magnitudes of
  218. // the components of the gradient is large).
  219. //
  220. // The original origin of this rescaling trick is somewhat unclear, the
  221. // earliest reference appears to be Oren [1], however it is widely
  222. // discussed without specific attribution in various texts including
  223. // [2] (p143).
  224. //
  225. // [1] Oren S.S., Self-scaling variable metric (SSVM) algorithms
  226. // Part II: Implementation and experiments, Management Science,
  227. // 20(5), 863-874, 1974.
  228. // [2] Nocedal J., Wright S., Numerical Optimization, Springer, 1999.
  229. const double approximate_eigenvalue_scale =
  230. delta_x_dot_delta_gradient / delta_gradient.dot(delta_gradient);
  231. inverse_hessian_ *= approximate_eigenvalue_scale;
  232. VLOG(4) << "Applying approximate_eigenvalue_scale: "
  233. << approximate_eigenvalue_scale << " to initial inverse "
  234. << "Hessian approximation.";
  235. }
  236. initialized_ = true;
  237. // Efficient O(num_parameters^2) BFGS update [2].
  238. //
  239. // Starting from dense BFGS update detailed in Nocedal [2] p140/177 and
  240. // using: y_k = delta_gradient, s_k = delta_x:
  241. //
  242. // \rho_k = 1.0 / (s_k' * y_k)
  243. // V_k = I - \rho_k * y_k * s_k'
  244. // H_k = (V_k' * H_{k-1} * V_k) + (\rho_k * s_k * s_k')
  245. //
  246. // This update involves matrix, matrix products which naively O(N^3),
  247. // however we can exploit our knowledge that H_k is positive definite
  248. // and thus by defn. symmetric to reduce the cost of the update:
  249. //
  250. // Expanding the update above yields:
  251. //
  252. // H_k = H_{k-1} +
  253. // \rho_k * ( (1.0 + \rho_k * y_k' * H_k * y_k) * s_k * s_k' -
  254. // (s_k * y_k' * H_k + H_k * y_k * s_k') )
  255. //
  256. // Using: A = (s_k * y_k' * H_k), and the knowledge that H_k = H_k', the
  257. // last term simplifies to (A + A'). Note that although A is not symmetric
  258. // (A + A') is symmetric. For ease of construction we also define
  259. // B = (1 + \rho_k * y_k' * H_k * y_k) * s_k * s_k', which is by defn
  260. // symmetric due to construction from: s_k * s_k'.
  261. //
  262. // Now we can write the BFGS update as:
  263. //
  264. // H_k = H_{k-1} + \rho_k * (B - (A + A'))
  265. // For efficiency, as H_k is by defn. symmetric, we will only maintain the
  266. // *lower* triangle of H_k (and all intermediary terms).
  267. const double rho_k = 1.0 / delta_x_dot_delta_gradient;
  268. // Calculate: A = s_k * y_k' * H_k
  269. Matrix A = delta_x * (delta_gradient.transpose() *
  270. inverse_hessian_.selfadjointView<Eigen::Lower>());
  271. // Calculate scalar: (1 + \rho_k * y_k' * H_k * y_k)
  272. const double delta_x_times_delta_x_transpose_scale_factor =
  273. (1.0 +
  274. (rho_k * delta_gradient.transpose() *
  275. inverse_hessian_.selfadjointView<Eigen::Lower>() * delta_gradient));
  276. // Calculate: B = (1 + \rho_k * y_k' * H_k * y_k) * s_k * s_k'
  277. Matrix B = Matrix::Zero(num_parameters_, num_parameters_);
  278. B.selfadjointView<Eigen::Lower>().rankUpdate(
  279. delta_x, delta_x_times_delta_x_transpose_scale_factor);
  280. // Finally, update inverse Hessian approximation according to:
  281. // H_k = H_{k-1} + \rho_k * (B - (A + A')). Note that (A + A') is
  282. // symmetric, even though A is not.
  283. inverse_hessian_.triangularView<Eigen::Lower>() +=
  284. rho_k * (B - A - A.transpose());
  285. }
  286. *search_direction = inverse_hessian_.selfadjointView<Eigen::Lower>() *
  287. (-1.0 * current.gradient);
  288. if (search_direction->dot(current.gradient) >= 0.0) {
  289. LOG(WARNING) << "Numerical failure in BFGS update: inverse Hessian "
  290. << "approximation is not positive definite, and thus "
  291. << "initial gradient for search direction is positive: "
  292. << search_direction->dot(current.gradient);
  293. is_positive_definite_ = false;
  294. return false;
  295. }
  296. return true;
  297. }
  298. private:
  299. const int num_parameters_;
  300. const bool use_approximate_eigenvalue_scaling_;
  301. Matrix inverse_hessian_;
  302. bool initialized_;
  303. bool is_positive_definite_;
  304. };
  305. LineSearchDirection::~LineSearchDirection() = default;
  306. std::unique_ptr<LineSearchDirection> LineSearchDirection::Create(
  307. const LineSearchDirection::Options& options) {
  308. if (options.type == STEEPEST_DESCENT) {
  309. return std::make_unique<SteepestDescent>();
  310. }
  311. if (options.type == NONLINEAR_CONJUGATE_GRADIENT) {
  312. return std::make_unique<NonlinearConjugateGradient>(
  313. options.nonlinear_conjugate_gradient_type, options.function_tolerance);
  314. }
  315. if (options.type == ceres::LBFGS) {
  316. return std::make_unique<ceres::internal::LBFGS>(
  317. options.num_parameters,
  318. options.max_lbfgs_rank,
  319. options.use_approximate_eigenvalue_bfgs_scaling);
  320. }
  321. if (options.type == ceres::BFGS) {
  322. return std::make_unique<ceres::internal::BFGS>(
  323. options.num_parameters,
  324. options.use_approximate_eigenvalue_bfgs_scaling);
  325. }
  326. LOG(ERROR) << "Unknown line search direction type: " << options.type;
  327. return nullptr;
  328. }
  329. } // namespace ceres::internal