low_rank_inverse_hessian.cc 7.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179
  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/low_rank_inverse_hessian.h"
  31. #include <list>
  32. #include "ceres/internal/eigen.h"
  33. #include "glog/logging.h"
  34. namespace ceres::internal {
  35. // The (L)BFGS algorithm explicitly requires that the secant equation:
  36. //
  37. // B_{k+1} * s_k = y_k
  38. //
  39. // Is satisfied at each iteration, where B_{k+1} is the approximated
  40. // Hessian at the k+1-th iteration, s_k = (x_{k+1} - x_{k}) and
  41. // y_k = (grad_{k+1} - grad_{k}). As the approximated Hessian must be
  42. // positive definite, this is equivalent to the condition:
  43. //
  44. // s_k^T * y_k > 0 [s_k^T * B_{k+1} * s_k = s_k^T * y_k > 0]
  45. //
  46. // This condition would always be satisfied if the function was strictly
  47. // convex, alternatively, it is always satisfied provided that a Wolfe line
  48. // search is used (even if the function is not strictly convex). See [1]
  49. // (p138) for a proof.
  50. //
  51. // Although Ceres will always use a Wolfe line search when using (L)BFGS,
  52. // practical implementation considerations mean that the line search
  53. // may return a point that satisfies only the Armijo condition, and thus
  54. // could violate the Secant equation. As such, we will only use a step
  55. // to update the Hessian approximation if:
  56. //
  57. // s_k^T * y_k > tolerance
  58. //
  59. // It is important that tolerance is very small (and >=0), as otherwise we
  60. // might skip the update too often and fail to capture important curvature
  61. // information in the Hessian. For example going from 1e-10 -> 1e-14 improves
  62. // the NIST benchmark score from 43/54 to 53/54.
  63. //
  64. // [1] Nocedal J., Wright S., Numerical Optimization, 2nd Ed. Springer, 1999.
  65. //
  66. // TODO(alexs.mac): Consider using Damped BFGS update instead of
  67. // skipping update.
  68. const double kLBFGSSecantConditionHessianUpdateTolerance = 1e-14;
  69. LowRankInverseHessian::LowRankInverseHessian(
  70. int num_parameters,
  71. int max_num_corrections,
  72. bool use_approximate_eigenvalue_scaling)
  73. : num_parameters_(num_parameters),
  74. max_num_corrections_(max_num_corrections),
  75. use_approximate_eigenvalue_scaling_(use_approximate_eigenvalue_scaling),
  76. approximate_eigenvalue_scale_(1.0),
  77. delta_x_history_(num_parameters, max_num_corrections),
  78. delta_gradient_history_(num_parameters, max_num_corrections),
  79. delta_x_dot_delta_gradient_(max_num_corrections) {}
  80. bool LowRankInverseHessian::Update(const Vector& delta_x,
  81. const Vector& delta_gradient) {
  82. const double delta_x_dot_delta_gradient = delta_x.dot(delta_gradient);
  83. if (delta_x_dot_delta_gradient <=
  84. kLBFGSSecantConditionHessianUpdateTolerance) {
  85. VLOG(2) << "Skipping L-BFGS Update, delta_x_dot_delta_gradient too "
  86. << "small: " << delta_x_dot_delta_gradient
  87. << ", tolerance: " << kLBFGSSecantConditionHessianUpdateTolerance
  88. << " (Secant condition).";
  89. return false;
  90. }
  91. int next = indices_.size();
  92. // Once the size of the list reaches max_num_corrections_, simulate
  93. // a circular buffer by removing the first element of the list and
  94. // making it the next position where the LBFGS history is stored.
  95. if (next == max_num_corrections_) {
  96. next = indices_.front();
  97. indices_.pop_front();
  98. }
  99. indices_.push_back(next);
  100. delta_x_history_.col(next) = delta_x;
  101. delta_gradient_history_.col(next) = delta_gradient;
  102. delta_x_dot_delta_gradient_(next) = delta_x_dot_delta_gradient;
  103. approximate_eigenvalue_scale_ =
  104. delta_x_dot_delta_gradient / delta_gradient.squaredNorm();
  105. return true;
  106. }
  107. void LowRankInverseHessian::RightMultiplyAndAccumulate(const double* x_ptr,
  108. double* y_ptr) const {
  109. ConstVectorRef gradient(x_ptr, num_parameters_);
  110. VectorRef search_direction(y_ptr, num_parameters_);
  111. search_direction = gradient;
  112. const int num_corrections = indices_.size();
  113. Vector alpha(num_corrections);
  114. for (auto it = indices_.rbegin(); it != indices_.rend(); ++it) {
  115. const double alpha_i = delta_x_history_.col(*it).dot(search_direction) /
  116. delta_x_dot_delta_gradient_(*it);
  117. search_direction -= alpha_i * delta_gradient_history_.col(*it);
  118. alpha(*it) = alpha_i;
  119. }
  120. if (use_approximate_eigenvalue_scaling_) {
  121. // Rescale the initial inverse Hessian approximation (H_0) to be iteratively
  122. // updated so that it is of similar 'size' to the true inverse Hessian along
  123. // the most recent search direction. As shown in [1]:
  124. //
  125. // \gamma_k = (delta_gradient_{k-1}' * delta_x_{k-1}) /
  126. // (delta_gradient_{k-1}' * delta_gradient_{k-1})
  127. //
  128. // Satisfies:
  129. //
  130. // (1 / \lambda_m) <= \gamma_k <= (1 / \lambda_1)
  131. //
  132. // Where \lambda_1 & \lambda_m are the smallest and largest eigenvalues of
  133. // the true Hessian (not the inverse) along the most recent search direction
  134. // respectively. Thus \gamma is an approximate eigenvalue of the true
  135. // inverse Hessian, and choosing: H_0 = I * \gamma will yield a starting
  136. // point that has a similar scale to the true inverse Hessian. This
  137. // technique is widely reported to often improve convergence, however this
  138. // is not universally true, particularly if there are errors in the initial
  139. // jacobians, or if there are significant differences in the sensitivity
  140. // of the problem to the parameters (i.e. the range of the magnitudes of
  141. // the components of the gradient is large).
  142. //
  143. // The original origin of this rescaling trick is somewhat unclear, the
  144. // earliest reference appears to be Oren [1], however it is widely discussed
  145. // without specific attribution in various texts including [2] (p143/178).
  146. //
  147. // [1] Oren S.S., Self-scaling variable metric (SSVM) algorithms Part II:
  148. // Implementation and experiments, Management Science,
  149. // 20(5), 863-874, 1974.
  150. // [2] Nocedal J., Wright S., Numerical Optimization, Springer, 1999.
  151. search_direction *= approximate_eigenvalue_scale_;
  152. VLOG(4) << "Applying approximate_eigenvalue_scale: "
  153. << approximate_eigenvalue_scale_ << " to initial inverse Hessian "
  154. << "approximation.";
  155. }
  156. for (const int i : indices_) {
  157. const double beta = delta_gradient_history_.col(i).dot(search_direction) /
  158. delta_x_dot_delta_gradient_(i);
  159. search_direction += delta_x_history_.col(i) * (alpha(i) - beta);
  160. }
  161. }
  162. } // namespace ceres::internal