conjugate_gradients_solver.h 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2022 Google Inc. All rights reserved.
  3. // http://ceres-solver.org/
  4. //
  5. // Redistribution and use in source and binary forms, with or without
  6. // modification, are permitted provided that the following conditions are met:
  7. //
  8. // * Redistributions of source code must retain the above copyright notice,
  9. // this list of conditions and the following disclaimer.
  10. // * Redistributions in binary form must reproduce the above copyright notice,
  11. // this list of conditions and the following disclaimer in the documentation
  12. // and/or other materials provided with the distribution.
  13. // * Neither the name of Google Inc. nor the names of its contributors may be
  14. // used to endorse or promote products derived from this software without
  15. // specific prior written permission.
  16. //
  17. // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
  18. // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
  19. // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
  20. // ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
  21. // LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
  22. // CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
  23. // SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
  24. // INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
  25. // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
  26. // ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
  27. // POSSIBILITY OF SUCH DAMAGE.
  28. //
  29. // Author: sameeragarwal@google.com (Sameer Agarwal)
  30. //
  31. // Preconditioned Conjugate Gradients based solver for positive
  32. // semidefinite linear systems.
  33. #ifndef CERES_INTERNAL_CONJUGATE_GRADIENTS_SOLVER_H_
  34. #define CERES_INTERNAL_CONJUGATE_GRADIENTS_SOLVER_H_
  35. #include <cmath>
  36. #include <cstddef>
  37. #include <utility>
  38. #include "ceres/eigen_vector_ops.h"
  39. #include "ceres/internal/disable_warnings.h"
  40. #include "ceres/internal/eigen.h"
  41. #include "ceres/internal/export.h"
  42. #include "ceres/linear_operator.h"
  43. #include "ceres/linear_solver.h"
  44. #include "ceres/stringprintf.h"
  45. #include "ceres/types.h"
  46. #include "glog/logging.h"
  47. namespace ceres::internal {
  48. // Interface for the linear operator used by ConjugateGradientsSolver.
  49. template <typename DenseVectorType>
  50. class ConjugateGradientsLinearOperator {
  51. public:
  52. ~ConjugateGradientsLinearOperator() = default;
  53. virtual void RightMultiplyAndAccumulate(const DenseVectorType& x,
  54. DenseVectorType& y) = 0;
  55. };
  56. // Adapter class that makes LinearOperator appear like an instance of
  57. // ConjugateGradientsLinearOperator.
  58. class LinearOperatorAdapter : public ConjugateGradientsLinearOperator<Vector> {
  59. public:
  60. LinearOperatorAdapter(LinearOperator& linear_operator)
  61. : linear_operator_(linear_operator) {}
  62. void RightMultiplyAndAccumulate(const Vector& x, Vector& y) final {
  63. linear_operator_.RightMultiplyAndAccumulate(x, y);
  64. }
  65. private:
  66. LinearOperator& linear_operator_;
  67. };
  68. // Options to control the ConjugateGradientsSolver. For detailed documentation
  69. // for each of these options see linear_solver.h
  70. struct ConjugateGradientsSolverOptions {
  71. int min_num_iterations = 1;
  72. int max_num_iterations = 1;
  73. int residual_reset_period = 10;
  74. double r_tolerance = 0.0;
  75. double q_tolerance = 0.0;
  76. ContextImpl* context = nullptr;
  77. int num_threads = 1;
  78. };
  79. // This function implements the now classical Conjugate Gradients algorithm of
  80. // Hestenes & Stiefel for solving positive semidefinite linear systems.
  81. // Optionally it can use a preconditioner also to reduce the condition number of
  82. // the linear system and improve the convergence rate. Modern references for
  83. // Conjugate Gradients are the books by Yousef Saad and Trefethen & Bau. This
  84. // implementation of CG has been augmented with additional termination tests
  85. // that are needed for forcing early termination when used as part of an inexact
  86. // Newton solver.
  87. //
  88. // This implementation is templated over DenseVectorType and then in turn on
  89. // ConjugateGradientsLinearOperator, which allows us to write an abstract
  90. // implementaion of the Conjugate Gradients algorithm without worrying about how
  91. // these objects are implemented or where they are stored. In particular it
  92. // allows us to have a single implementation that works on CPU and GPU based
  93. // matrices and vectors.
  94. //
  95. // scratch must contain pointers to four DenseVector objects of the same size as
  96. // rhs and solution. By asking the user for scratch space, we guarantee that we
  97. // will not perform any allocations inside this function.
  98. template <typename DenseVectorType>
  99. LinearSolver::Summary ConjugateGradientsSolver(
  100. const ConjugateGradientsSolverOptions options,
  101. ConjugateGradientsLinearOperator<DenseVectorType>& lhs,
  102. const DenseVectorType& rhs,
  103. ConjugateGradientsLinearOperator<DenseVectorType>& preconditioner,
  104. DenseVectorType* scratch[4],
  105. DenseVectorType& solution) {
  106. auto IsZeroOrInfinity = [](double x) {
  107. return ((x == 0.0) || std::isinf(x));
  108. };
  109. DenseVectorType& p = *scratch[0];
  110. DenseVectorType& r = *scratch[1];
  111. DenseVectorType& z = *scratch[2];
  112. DenseVectorType& tmp = *scratch[3];
  113. LinearSolver::Summary summary;
  114. summary.termination_type = LinearSolverTerminationType::NO_CONVERGENCE;
  115. summary.message = "Maximum number of iterations reached.";
  116. summary.num_iterations = 0;
  117. const double norm_rhs = Norm(rhs, options.context, options.num_threads);
  118. if (norm_rhs == 0.0) {
  119. SetZero(solution, options.context, options.num_threads);
  120. summary.termination_type = LinearSolverTerminationType::SUCCESS;
  121. summary.message = "Convergence. |b| = 0.";
  122. return summary;
  123. }
  124. const double tol_r = options.r_tolerance * norm_rhs;
  125. SetZero(tmp, options.context, options.num_threads);
  126. lhs.RightMultiplyAndAccumulate(solution, tmp);
  127. // r = rhs - tmp
  128. Axpby(1.0, rhs, -1.0, tmp, r, options.context, options.num_threads);
  129. double norm_r = Norm(r, options.context, options.num_threads);
  130. if (options.min_num_iterations == 0 && norm_r <= tol_r) {
  131. summary.termination_type = LinearSolverTerminationType::SUCCESS;
  132. summary.message =
  133. StringPrintf("Convergence. |r| = %e <= %e.", norm_r, tol_r);
  134. return summary;
  135. }
  136. double rho = 1.0;
  137. // Initial value of the quadratic model Q = x'Ax - 2 * b'x.
  138. // double Q0 = -1.0 * solution.dot(rhs + r);
  139. Axpby(1.0, rhs, 1.0, r, tmp, options.context, options.num_threads);
  140. double Q0 = -Dot(solution, tmp, options.context, options.num_threads);
  141. for (summary.num_iterations = 1;; ++summary.num_iterations) {
  142. SetZero(z, options.context, options.num_threads);
  143. preconditioner.RightMultiplyAndAccumulate(r, z);
  144. const double last_rho = rho;
  145. // rho = r.dot(z);
  146. rho = Dot(r, z, options.context, options.num_threads);
  147. if (IsZeroOrInfinity(rho)) {
  148. summary.termination_type = LinearSolverTerminationType::FAILURE;
  149. summary.message = StringPrintf("Numerical failure. rho = r'z = %e.", rho);
  150. break;
  151. }
  152. if (summary.num_iterations == 1) {
  153. Copy(z, p, options.context, options.num_threads);
  154. } else {
  155. const double beta = rho / last_rho;
  156. if (IsZeroOrInfinity(beta)) {
  157. summary.termination_type = LinearSolverTerminationType::FAILURE;
  158. summary.message = StringPrintf(
  159. "Numerical failure. beta = rho_n / rho_{n-1} = %e, "
  160. "rho_n = %e, rho_{n-1} = %e",
  161. beta,
  162. rho,
  163. last_rho);
  164. break;
  165. }
  166. // p = z + beta * p;
  167. Axpby(1.0, z, beta, p, p, options.context, options.num_threads);
  168. }
  169. DenseVectorType& q = z;
  170. SetZero(q, options.context, options.num_threads);
  171. lhs.RightMultiplyAndAccumulate(p, q);
  172. const double pq = Dot(p, q, options.context, options.num_threads);
  173. if ((pq <= 0) || std::isinf(pq)) {
  174. summary.termination_type = LinearSolverTerminationType::NO_CONVERGENCE;
  175. summary.message = StringPrintf(
  176. "Matrix is indefinite, no more progress can be made. "
  177. "p'q = %e. |p| = %e, |q| = %e",
  178. pq,
  179. Norm(p, options.context, options.num_threads),
  180. Norm(q, options.context, options.num_threads));
  181. break;
  182. }
  183. const double alpha = rho / pq;
  184. if (std::isinf(alpha)) {
  185. summary.termination_type = LinearSolverTerminationType::FAILURE;
  186. summary.message = StringPrintf(
  187. "Numerical failure. alpha = rho / pq = %e, rho = %e, pq = %e.",
  188. alpha,
  189. rho,
  190. pq);
  191. break;
  192. }
  193. // solution = solution + alpha * p;
  194. Axpby(1.0,
  195. solution,
  196. alpha,
  197. p,
  198. solution,
  199. options.context,
  200. options.num_threads);
  201. // Ideally we would just use the update r = r - alpha*q to keep
  202. // track of the residual vector. However this estimate tends to
  203. // drift over time due to round off errors. Thus every
  204. // residual_reset_period iterations, we calculate the residual as
  205. // r = b - Ax. We do not do this every iteration because this
  206. // requires an additional matrix vector multiply which would
  207. // double the complexity of the CG algorithm.
  208. if (summary.num_iterations % options.residual_reset_period == 0) {
  209. SetZero(tmp, options.context, options.num_threads);
  210. lhs.RightMultiplyAndAccumulate(solution, tmp);
  211. Axpby(1.0, rhs, -1.0, tmp, r, options.context, options.num_threads);
  212. // r = rhs - tmp;
  213. } else {
  214. Axpby(1.0, r, -alpha, q, r, options.context, options.num_threads);
  215. // r = r - alpha * q;
  216. }
  217. // Quadratic model based termination.
  218. // Q1 = x'Ax - 2 * b' x.
  219. // const double Q1 = -1.0 * solution.dot(rhs + r);
  220. Axpby(1.0, rhs, 1.0, r, tmp, options.context, options.num_threads);
  221. const double Q1 = -Dot(solution, tmp, options.context, options.num_threads);
  222. // For PSD matrices A, let
  223. //
  224. // Q(x) = x'Ax - 2b'x
  225. //
  226. // be the cost of the quadratic function defined by A and b. Then,
  227. // the solver terminates at iteration i if
  228. //
  229. // i * (Q(x_i) - Q(x_i-1)) / Q(x_i) < q_tolerance.
  230. //
  231. // This termination criterion is more useful when using CG to
  232. // solve the Newton step. This particular convergence test comes
  233. // from Stephen Nash's work on truncated Newton
  234. // methods. References:
  235. //
  236. // 1. Stephen G. Nash & Ariela Sofer, Assessing A Search
  237. // Direction Within A Truncated Newton Method, Operation
  238. // Research Letters 9(1990) 219-221.
  239. //
  240. // 2. Stephen G. Nash, A Survey of Truncated Newton Methods,
  241. // Journal of Computational and Applied Mathematics,
  242. // 124(1-2), 45-59, 2000.
  243. //
  244. const double zeta = summary.num_iterations * (Q1 - Q0) / Q1;
  245. if (zeta < options.q_tolerance &&
  246. summary.num_iterations >= options.min_num_iterations) {
  247. summary.termination_type = LinearSolverTerminationType::SUCCESS;
  248. summary.message =
  249. StringPrintf("Iteration: %d Convergence: zeta = %e < %e. |r| = %e",
  250. summary.num_iterations,
  251. zeta,
  252. options.q_tolerance,
  253. Norm(r, options.context, options.num_threads));
  254. break;
  255. }
  256. Q0 = Q1;
  257. // Residual based termination.
  258. norm_r = Norm(r, options.context, options.num_threads);
  259. if (norm_r <= tol_r &&
  260. summary.num_iterations >= options.min_num_iterations) {
  261. summary.termination_type = LinearSolverTerminationType::SUCCESS;
  262. summary.message =
  263. StringPrintf("Iteration: %d Convergence. |r| = %e <= %e.",
  264. summary.num_iterations,
  265. norm_r,
  266. tol_r);
  267. break;
  268. }
  269. if (summary.num_iterations >= options.max_num_iterations) {
  270. break;
  271. }
  272. }
  273. return summary;
  274. }
  275. } // namespace ceres::internal
  276. #include "ceres/internal/reenable_warnings.h"
  277. #endif // CERES_INTERNAL_CONJUGATE_GRADIENTS_SOLVER_H_