schur_complement_solver.h 7.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2015 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. #ifndef CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_
  31. #define CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_
  32. #include <memory>
  33. #include <set>
  34. #include <utility>
  35. #include <vector>
  36. #include "ceres/block_random_access_diagonal_matrix.h"
  37. #include "ceres/block_random_access_matrix.h"
  38. #include "ceres/block_sparse_matrix.h"
  39. #include "ceres/block_structure.h"
  40. #include "ceres/dense_cholesky.h"
  41. #include "ceres/internal/config.h"
  42. #include "ceres/internal/export.h"
  43. #include "ceres/linear_solver.h"
  44. #include "ceres/schur_eliminator.h"
  45. #include "ceres/types.h"
  46. #ifdef CERES_USE_EIGEN_SPARSE
  47. #include "Eigen/OrderingMethods"
  48. #include "Eigen/SparseCholesky"
  49. #endif
  50. #include "ceres/internal/disable_warnings.h"
  51. namespace ceres::internal {
  52. class BlockSparseMatrix;
  53. class SparseCholesky;
  54. // Base class for Schur complement based linear least squares
  55. // solvers. It assumes that the input linear system Ax = b can be
  56. // partitioned into
  57. //
  58. // E y + F z = b
  59. //
  60. // Where x = [y;z] is a partition of the variables. The partitioning
  61. // of the variables is such that, E'E is a block diagonal
  62. // matrix. Further, the rows of A are ordered so that for every
  63. // variable block in y, all the rows containing that variable block
  64. // occur as a vertically contiguous block. i.e the matrix A looks like
  65. //
  66. // E F
  67. // A = [ y1 0 0 0 | z1 0 0 0 z5]
  68. // [ y1 0 0 0 | z1 z2 0 0 0]
  69. // [ 0 y2 0 0 | 0 0 z3 0 0]
  70. // [ 0 0 y3 0 | z1 z2 z3 z4 z5]
  71. // [ 0 0 y3 0 | z1 0 0 0 z5]
  72. // [ 0 0 0 y4 | 0 0 0 0 z5]
  73. // [ 0 0 0 y4 | 0 z2 0 0 0]
  74. // [ 0 0 0 y4 | 0 0 0 0 0]
  75. // [ 0 0 0 0 | z1 0 0 0 0]
  76. // [ 0 0 0 0 | 0 0 z3 z4 z5]
  77. //
  78. // This structure should be reflected in the corresponding
  79. // CompressedRowBlockStructure object associated with A. The linear
  80. // system Ax = b should either be well posed or the array D below
  81. // should be non-null and the diagonal matrix corresponding to it
  82. // should be non-singular.
  83. //
  84. // SchurComplementSolver has two sub-classes.
  85. //
  86. // DenseSchurComplementSolver: For problems where the Schur complement
  87. // matrix is small and dense, or if CHOLMOD/SuiteSparse is not
  88. // installed. For structure from motion problems, this is solver can
  89. // be used for problems with upto a few hundred cameras.
  90. //
  91. // SparseSchurComplementSolver: For problems where the Schur
  92. // complement matrix is large and sparse. It requires that Ceres be
  93. // build with at least one sparse linear algebra library, as it
  94. // computes a sparse Cholesky factorization of the Schur complement.
  95. //
  96. // This solver can be used for solving structure from motion problems
  97. // with tens of thousands of cameras, though depending on the exact
  98. // sparsity structure, it maybe better to use an iterative solver.
  99. //
  100. // The two solvers can be instantiated by calling
  101. // LinearSolver::CreateLinearSolver with LinearSolver::Options::type
  102. // set to DENSE_SCHUR and SPARSE_SCHUR
  103. // respectively. LinearSolver::Options::elimination_groups[0] should
  104. // be at least 1.
  105. class CERES_NO_EXPORT SchurComplementSolver : public BlockSparseMatrixSolver {
  106. public:
  107. explicit SchurComplementSolver(const LinearSolver::Options& options);
  108. SchurComplementSolver(const SchurComplementSolver&) = delete;
  109. void operator=(const SchurComplementSolver&) = delete;
  110. LinearSolver::Summary SolveImpl(
  111. BlockSparseMatrix* A,
  112. const double* b,
  113. const LinearSolver::PerSolveOptions& per_solve_options,
  114. double* x) override;
  115. protected:
  116. const LinearSolver::Options& options() const { return options_; }
  117. void set_lhs(std::unique_ptr<BlockRandomAccessMatrix> lhs) {
  118. lhs_ = std::move(lhs);
  119. }
  120. const BlockRandomAccessMatrix* lhs() const { return lhs_.get(); }
  121. BlockRandomAccessMatrix* mutable_lhs() { return lhs_.get(); }
  122. void ResizeRhs(int n) { rhs_.resize(n); }
  123. const Vector& rhs() const { return rhs_; }
  124. private:
  125. virtual void InitStorage(const CompressedRowBlockStructure* bs) = 0;
  126. virtual LinearSolver::Summary SolveReducedLinearSystem(
  127. const LinearSolver::PerSolveOptions& per_solve_options,
  128. double* solution) = 0;
  129. LinearSolver::Options options_;
  130. std::unique_ptr<SchurEliminatorBase> eliminator_;
  131. std::unique_ptr<BlockRandomAccessMatrix> lhs_;
  132. Vector rhs_;
  133. };
  134. // Dense Cholesky factorization based solver.
  135. class CERES_NO_EXPORT DenseSchurComplementSolver final
  136. : public SchurComplementSolver {
  137. public:
  138. explicit DenseSchurComplementSolver(const LinearSolver::Options& options);
  139. DenseSchurComplementSolver(const DenseSchurComplementSolver&) = delete;
  140. void operator=(const DenseSchurComplementSolver&) = delete;
  141. ~DenseSchurComplementSolver() override;
  142. private:
  143. void InitStorage(const CompressedRowBlockStructure* bs) final;
  144. LinearSolver::Summary SolveReducedLinearSystem(
  145. const LinearSolver::PerSolveOptions& per_solve_options,
  146. double* solution) final;
  147. std::unique_ptr<DenseCholesky> cholesky_;
  148. };
  149. // Sparse Cholesky factorization based solver.
  150. class CERES_NO_EXPORT SparseSchurComplementSolver final
  151. : public SchurComplementSolver {
  152. public:
  153. explicit SparseSchurComplementSolver(const LinearSolver::Options& options);
  154. SparseSchurComplementSolver(const SparseSchurComplementSolver&) = delete;
  155. void operator=(const SparseSchurComplementSolver&) = delete;
  156. ~SparseSchurComplementSolver() override;
  157. private:
  158. void InitStorage(const CompressedRowBlockStructure* bs) final;
  159. LinearSolver::Summary SolveReducedLinearSystem(
  160. const LinearSolver::PerSolveOptions& per_solve_options,
  161. double* solution) final;
  162. LinearSolver::Summary SolveReducedLinearSystemUsingConjugateGradients(
  163. const LinearSolver::PerSolveOptions& per_solve_options, double* solution);
  164. std::vector<Block> blocks_;
  165. std::unique_ptr<SparseCholesky> sparse_cholesky_;
  166. std::unique_ptr<BlockRandomAccessDiagonalMatrix> preconditioner_;
  167. std::unique_ptr<CompressedRowSparseMatrix> crs_lhs_;
  168. Vector cg_solution_;
  169. Vector* scratch_[4] = {nullptr, nullptr, nullptr, nullptr};
  170. };
  171. } // namespace ceres::internal
  172. #include "ceres/internal/reenable_warnings.h"
  173. #endif // CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_