implicit_schur_complement_test.cc 8.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240
  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/implicit_schur_complement.h"
  31. #include <cstddef>
  32. #include <memory>
  33. #include "Eigen/Dense"
  34. #include "ceres/block_random_access_dense_matrix.h"
  35. #include "ceres/block_sparse_matrix.h"
  36. #include "ceres/casts.h"
  37. #include "ceres/context_impl.h"
  38. #include "ceres/internal/eigen.h"
  39. #include "ceres/linear_least_squares_problems.h"
  40. #include "ceres/linear_solver.h"
  41. #include "ceres/schur_eliminator.h"
  42. #include "ceres/triplet_sparse_matrix.h"
  43. #include "ceres/types.h"
  44. #include "glog/logging.h"
  45. #include "gtest/gtest.h"
  46. namespace ceres::internal {
  47. using testing::AssertionResult;
  48. const double kEpsilon = 1e-14;
  49. class ImplicitSchurComplementTest : public ::testing::Test {
  50. protected:
  51. void SetUp() final {
  52. auto problem = CreateLinearLeastSquaresProblemFromId(2);
  53. CHECK(problem != nullptr);
  54. A_.reset(down_cast<BlockSparseMatrix*>(problem->A.release()));
  55. b_ = std::move(problem->b);
  56. D_ = std::move(problem->D);
  57. num_cols_ = A_->num_cols();
  58. num_rows_ = A_->num_rows();
  59. num_eliminate_blocks_ = problem->num_eliminate_blocks;
  60. }
  61. void ReducedLinearSystemAndSolution(double* D,
  62. Matrix* lhs,
  63. Vector* rhs,
  64. Vector* solution) {
  65. const CompressedRowBlockStructure* bs = A_->block_structure();
  66. const int num_col_blocks = bs->cols.size();
  67. auto blocks = Tail(bs->cols, num_col_blocks - num_eliminate_blocks_);
  68. BlockRandomAccessDenseMatrix blhs(blocks, &context_, 1);
  69. const int num_schur_rows = blhs.num_rows();
  70. LinearSolver::Options options;
  71. options.elimination_groups.push_back(num_eliminate_blocks_);
  72. options.type = DENSE_SCHUR;
  73. ContextImpl context;
  74. options.context = &context;
  75. std::unique_ptr<SchurEliminatorBase> eliminator =
  76. SchurEliminatorBase::Create(options);
  77. CHECK(eliminator != nullptr);
  78. const bool kFullRankETE = true;
  79. eliminator->Init(num_eliminate_blocks_, kFullRankETE, bs);
  80. lhs->resize(num_schur_rows, num_schur_rows);
  81. rhs->resize(num_schur_rows);
  82. eliminator->Eliminate(
  83. BlockSparseMatrixData(*A_), b_.get(), D, &blhs, rhs->data());
  84. MatrixRef lhs_ref(blhs.mutable_values(), num_schur_rows, num_schur_rows);
  85. // lhs_ref is an upper triangular matrix. Construct a full version
  86. // of lhs_ref in lhs by transposing lhs_ref, choosing the strictly
  87. // lower triangular part of the matrix and adding it to lhs_ref.
  88. *lhs = lhs_ref;
  89. lhs->triangularView<Eigen::StrictlyLower>() =
  90. lhs_ref.triangularView<Eigen::StrictlyUpper>().transpose();
  91. solution->resize(num_cols_);
  92. solution->setZero();
  93. VectorRef schur_solution(solution->data() + num_cols_ - num_schur_rows,
  94. num_schur_rows);
  95. schur_solution = lhs->selfadjointView<Eigen::Upper>().llt().solve(*rhs);
  96. eliminator->BackSubstitute(BlockSparseMatrixData(*A_),
  97. b_.get(),
  98. D,
  99. schur_solution.data(),
  100. solution->data());
  101. }
  102. AssertionResult TestImplicitSchurComplement(double* D) {
  103. Matrix lhs;
  104. Vector rhs;
  105. Vector reference_solution;
  106. ReducedLinearSystemAndSolution(D, &lhs, &rhs, &reference_solution);
  107. LinearSolver::Options options;
  108. options.elimination_groups.push_back(num_eliminate_blocks_);
  109. options.preconditioner_type = JACOBI;
  110. ContextImpl context;
  111. options.context = &context;
  112. ImplicitSchurComplement isc(options);
  113. isc.Init(*A_, D, b_.get());
  114. const int num_f_cols = lhs.cols();
  115. const int num_e_cols = num_cols_ - num_f_cols;
  116. Matrix A_dense, E, F, DE, DF;
  117. A_->ToDenseMatrix(&A_dense);
  118. E = A_dense.leftCols(A_->num_cols() - num_f_cols);
  119. F = A_dense.rightCols(num_f_cols);
  120. if (D) {
  121. DE = VectorRef(D, num_e_cols).asDiagonal();
  122. DF = VectorRef(D + num_e_cols, num_f_cols).asDiagonal();
  123. } else {
  124. DE = Matrix::Zero(num_e_cols, num_e_cols);
  125. DF = Matrix::Zero(num_f_cols, num_f_cols);
  126. }
  127. // Z = (block_diagonal(F'F))^-1 F'E (E'E)^-1 E'F
  128. // Here, assuming that block_diagonal(F'F) == diagonal(F'F)
  129. Matrix Z_reference =
  130. (F.transpose() * F + DF).diagonal().asDiagonal().inverse() *
  131. F.transpose() * E * (E.transpose() * E + DE).inverse() * E.transpose() *
  132. F;
  133. for (int i = 0; i < num_f_cols; ++i) {
  134. Vector x(num_f_cols);
  135. x.setZero();
  136. x(i) = 1.0;
  137. Vector y(num_f_cols);
  138. y = lhs * x;
  139. Vector z(num_f_cols);
  140. isc.RightMultiplyAndAccumulate(x.data(), z.data());
  141. // The i^th column of the implicit schur complement is the same as
  142. // the explicit schur complement.
  143. if ((y - z).norm() > kEpsilon) {
  144. return testing::AssertionFailure()
  145. << "Explicit and Implicit SchurComplements differ in "
  146. << "column " << i << ". explicit: " << y.transpose()
  147. << " implicit: " << z.transpose();
  148. }
  149. y.setZero();
  150. y = Z_reference * x;
  151. z.setZero();
  152. isc.InversePowerSeriesOperatorRightMultiplyAccumulate(x.data(), z.data());
  153. // The i^th column of operator Z stored implicitly is the same as its
  154. // explicit version.
  155. if ((y - z).norm() > kEpsilon) {
  156. return testing::AssertionFailure()
  157. << "Explicit and Implicit operators used to approximate the "
  158. "inversion of schur complement via power series expansion "
  159. "differ in column "
  160. << i << ". explicit: " << y.transpose()
  161. << " implicit: " << z.transpose();
  162. }
  163. }
  164. // Compare the rhs of the reduced linear system
  165. if ((isc.rhs() - rhs).norm() > kEpsilon) {
  166. return testing::AssertionFailure()
  167. << "Explicit and Implicit SchurComplements differ in "
  168. << "rhs. explicit: " << rhs.transpose()
  169. << " implicit: " << isc.rhs().transpose();
  170. }
  171. // Reference solution to the f_block.
  172. const Vector reference_f_sol =
  173. lhs.selfadjointView<Eigen::Upper>().llt().solve(rhs);
  174. // Backsubstituted solution from the implicit schur solver using the
  175. // reference solution to the f_block.
  176. Vector sol(num_cols_);
  177. isc.BackSubstitute(reference_f_sol.data(), sol.data());
  178. if ((sol - reference_solution).norm() > kEpsilon) {
  179. return testing::AssertionFailure()
  180. << "Explicit and Implicit SchurComplements solutions differ. "
  181. << "explicit: " << reference_solution.transpose()
  182. << " implicit: " << sol.transpose();
  183. }
  184. return testing::AssertionSuccess();
  185. }
  186. ContextImpl context_;
  187. int num_rows_;
  188. int num_cols_;
  189. int num_eliminate_blocks_;
  190. std::unique_ptr<BlockSparseMatrix> A_;
  191. std::unique_ptr<double[]> b_;
  192. std::unique_ptr<double[]> D_;
  193. };
  194. // Verify that the Schur Complement matrix implied by the
  195. // ImplicitSchurComplement class matches the one explicitly computed
  196. // by the SchurComplement solver.
  197. //
  198. // We do this with and without regularization to check that the
  199. // support for the LM diagonal is correct.
  200. TEST_F(ImplicitSchurComplementTest, SchurMatrixValuesTest) {
  201. EXPECT_TRUE(TestImplicitSchurComplement(nullptr));
  202. EXPECT_TRUE(TestImplicitSchurComplement(D_.get()));
  203. }
  204. } // namespace ceres::internal