schur_eliminator_test.cc 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373
  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/schur_eliminator.h"
  31. #include <algorithm>
  32. #include <memory>
  33. #include <random>
  34. #include <vector>
  35. #include "Eigen/Dense"
  36. #include "ceres/block_random_access_dense_matrix.h"
  37. #include "ceres/block_sparse_matrix.h"
  38. #include "ceres/block_structure.h"
  39. #include "ceres/casts.h"
  40. #include "ceres/context_impl.h"
  41. #include "ceres/detect_structure.h"
  42. #include "ceres/internal/eigen.h"
  43. #include "ceres/linear_least_squares_problems.h"
  44. #include "ceres/test_util.h"
  45. #include "ceres/triplet_sparse_matrix.h"
  46. #include "ceres/types.h"
  47. #include "glog/logging.h"
  48. #include "gtest/gtest.h"
  49. // TODO(sameeragarwal): Reduce the size of these tests and redo the
  50. // parameterization to be more efficient.
  51. namespace ceres::internal {
  52. class SchurEliminatorTest : public ::testing::Test {
  53. protected:
  54. void SetUpFromId(int id) {
  55. auto problem = CreateLinearLeastSquaresProblemFromId(id);
  56. CHECK(problem != nullptr);
  57. SetupHelper(problem.get());
  58. }
  59. void SetupHelper(LinearLeastSquaresProblem* problem) {
  60. A.reset(down_cast<BlockSparseMatrix*>(problem->A.release()));
  61. b = std::move(problem->b);
  62. D = std::move(problem->D);
  63. num_eliminate_blocks = problem->num_eliminate_blocks;
  64. num_eliminate_cols = 0;
  65. const CompressedRowBlockStructure* bs = A->block_structure();
  66. for (int i = 0; i < num_eliminate_blocks; ++i) {
  67. num_eliminate_cols += bs->cols[i].size;
  68. }
  69. }
  70. // Compute the golden values for the reduced linear system and the
  71. // solution to the linear least squares problem using dense linear
  72. // algebra.
  73. void ComputeReferenceSolution(const Vector& D) {
  74. Matrix J;
  75. A->ToDenseMatrix(&J);
  76. VectorRef f(b.get(), J.rows());
  77. Matrix H = (D.cwiseProduct(D)).asDiagonal();
  78. H.noalias() += J.transpose() * J;
  79. const Vector g = J.transpose() * f;
  80. const int schur_size = J.cols() - num_eliminate_cols;
  81. lhs_expected.resize(schur_size, schur_size);
  82. lhs_expected.setZero();
  83. rhs_expected.resize(schur_size);
  84. rhs_expected.setZero();
  85. sol_expected.resize(J.cols());
  86. sol_expected.setZero();
  87. Matrix P = H.block(0, 0, num_eliminate_cols, num_eliminate_cols);
  88. Matrix Q = H.block(0, num_eliminate_cols, num_eliminate_cols, schur_size);
  89. Matrix R =
  90. H.block(num_eliminate_cols, num_eliminate_cols, schur_size, schur_size);
  91. int row = 0;
  92. const CompressedRowBlockStructure* bs = A->block_structure();
  93. for (int i = 0; i < num_eliminate_blocks; ++i) {
  94. const int block_size = bs->cols[i].size;
  95. P.block(row, row, block_size, block_size) =
  96. P.block(row, row, block_size, block_size)
  97. .llt()
  98. .solve(Matrix::Identity(block_size, block_size));
  99. row += block_size;
  100. }
  101. lhs_expected.triangularView<Eigen::Upper>() = R - Q.transpose() * P * Q;
  102. rhs_expected =
  103. g.tail(schur_size) - Q.transpose() * P * g.head(num_eliminate_cols);
  104. sol_expected = H.llt().solve(g);
  105. }
  106. void EliminateSolveAndCompare(const VectorRef& diagonal,
  107. bool use_static_structure,
  108. const double relative_tolerance) {
  109. const CompressedRowBlockStructure* bs = A->block_structure();
  110. const int num_col_blocks = bs->cols.size();
  111. auto blocks = Tail(bs->cols, num_col_blocks - num_eliminate_blocks);
  112. BlockRandomAccessDenseMatrix lhs(blocks, &context_, 1);
  113. const int num_cols = A->num_cols();
  114. const int schur_size = lhs.num_rows();
  115. Vector rhs(schur_size);
  116. LinearSolver::Options options;
  117. options.context = &context_;
  118. options.elimination_groups.push_back(num_eliminate_blocks);
  119. if (use_static_structure) {
  120. DetectStructure(*bs,
  121. num_eliminate_blocks,
  122. &options.row_block_size,
  123. &options.e_block_size,
  124. &options.f_block_size);
  125. }
  126. std::unique_ptr<SchurEliminatorBase> eliminator =
  127. SchurEliminatorBase::Create(options);
  128. const bool kFullRankETE = true;
  129. eliminator->Init(num_eliminate_blocks, kFullRankETE, A->block_structure());
  130. eliminator->Eliminate(
  131. BlockSparseMatrixData(*A), b.get(), diagonal.data(), &lhs, rhs.data());
  132. MatrixRef lhs_ref(lhs.mutable_values(), lhs.num_rows(), lhs.num_cols());
  133. Vector reduced_sol =
  134. lhs_ref.selfadjointView<Eigen::Upper>().llt().solve(rhs);
  135. // Solution to the linear least squares problem.
  136. Vector sol(num_cols);
  137. sol.setZero();
  138. sol.tail(schur_size) = reduced_sol;
  139. eliminator->BackSubstitute(BlockSparseMatrixData(*A),
  140. b.get(),
  141. diagonal.data(),
  142. reduced_sol.data(),
  143. sol.data());
  144. Matrix delta = (lhs_ref - lhs_expected).selfadjointView<Eigen::Upper>();
  145. double diff = delta.norm();
  146. EXPECT_NEAR(diff / lhs_expected.norm(), 0.0, relative_tolerance);
  147. EXPECT_NEAR((rhs - rhs_expected).norm() / rhs_expected.norm(),
  148. 0.0,
  149. relative_tolerance);
  150. EXPECT_NEAR((sol - sol_expected).norm() / sol_expected.norm(),
  151. 0.0,
  152. relative_tolerance);
  153. }
  154. ContextImpl context_;
  155. std::unique_ptr<BlockSparseMatrix> A;
  156. std::unique_ptr<double[]> b;
  157. std::unique_ptr<double[]> D;
  158. int num_eliminate_blocks;
  159. int num_eliminate_cols;
  160. Matrix lhs_expected;
  161. Vector rhs_expected;
  162. Vector sol_expected;
  163. };
  164. TEST_F(SchurEliminatorTest, ScalarProblemNoRegularization) {
  165. SetUpFromId(2);
  166. Vector zero(A->num_cols());
  167. zero.setZero();
  168. ComputeReferenceSolution(VectorRef(zero.data(), A->num_cols()));
  169. EliminateSolveAndCompare(VectorRef(zero.data(), A->num_cols()), true, 1e-14);
  170. EliminateSolveAndCompare(VectorRef(zero.data(), A->num_cols()), false, 1e-14);
  171. }
  172. TEST_F(SchurEliminatorTest, ScalarProblemWithRegularization) {
  173. SetUpFromId(2);
  174. ComputeReferenceSolution(VectorRef(D.get(), A->num_cols()));
  175. EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), true, 1e-14);
  176. EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), false, 1e-14);
  177. }
  178. TEST_F(SchurEliminatorTest, VaryingFBlockSizeWithStaticStructure) {
  179. SetUpFromId(4);
  180. ComputeReferenceSolution(VectorRef(D.get(), A->num_cols()));
  181. EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), true, 1e-14);
  182. }
  183. TEST_F(SchurEliminatorTest, VaryingFBlockSizeWithoutStaticStructure) {
  184. SetUpFromId(4);
  185. ComputeReferenceSolution(VectorRef(D.get(), A->num_cols()));
  186. EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), false, 1e-14);
  187. }
  188. TEST(SchurEliminatorForOneFBlock, MatchesSchurEliminator) {
  189. constexpr int kRowBlockSize = 2;
  190. constexpr int kEBlockSize = 3;
  191. constexpr int kFBlockSize = 6;
  192. constexpr int num_e_blocks = 5;
  193. ContextImpl context;
  194. auto* bs = new CompressedRowBlockStructure;
  195. bs->cols.resize(num_e_blocks + 1);
  196. int col_pos = 0;
  197. for (int i = 0; i < num_e_blocks; ++i) {
  198. bs->cols[i].position = col_pos;
  199. bs->cols[i].size = kEBlockSize;
  200. col_pos += kEBlockSize;
  201. }
  202. bs->cols.back().position = col_pos;
  203. bs->cols.back().size = kFBlockSize;
  204. bs->rows.resize(2 * num_e_blocks + 1);
  205. int row_pos = 0;
  206. int cell_pos = 0;
  207. for (int i = 0; i < num_e_blocks; ++i) {
  208. {
  209. auto& row = bs->rows[2 * i];
  210. row.block.position = row_pos;
  211. row.block.size = kRowBlockSize;
  212. row_pos += kRowBlockSize;
  213. auto& cells = row.cells;
  214. cells.resize(2);
  215. cells[0].block_id = i;
  216. cells[0].position = cell_pos;
  217. cell_pos += kRowBlockSize * kEBlockSize;
  218. cells[1].block_id = num_e_blocks;
  219. cells[1].position = cell_pos;
  220. cell_pos += kRowBlockSize * kFBlockSize;
  221. }
  222. {
  223. auto& row = bs->rows[2 * i + 1];
  224. row.block.position = row_pos;
  225. row.block.size = kRowBlockSize;
  226. row_pos += kRowBlockSize;
  227. auto& cells = row.cells;
  228. cells.resize(1);
  229. cells[0].block_id = i;
  230. cells[0].position = cell_pos;
  231. cell_pos += kRowBlockSize * kEBlockSize;
  232. }
  233. }
  234. {
  235. auto& row = bs->rows.back();
  236. row.block.position = row_pos;
  237. row.block.size = kEBlockSize;
  238. row_pos += kRowBlockSize;
  239. auto& cells = row.cells;
  240. cells.resize(1);
  241. cells[0].block_id = num_e_blocks;
  242. cells[0].position = cell_pos;
  243. cell_pos += kEBlockSize * kEBlockSize;
  244. }
  245. BlockSparseMatrix matrix(bs);
  246. double* values = matrix.mutable_values();
  247. std::mt19937 prng;
  248. std::normal_distribution<> standard_normal;
  249. std::generate_n(values, matrix.num_nonzeros(), [&prng, &standard_normal] {
  250. return standard_normal(prng);
  251. });
  252. Vector b(matrix.num_rows());
  253. b.setRandom();
  254. Vector diagonal(matrix.num_cols());
  255. diagonal.setOnes();
  256. std::vector<Block> blocks;
  257. blocks.emplace_back(kFBlockSize, 0);
  258. BlockRandomAccessDenseMatrix actual_lhs(blocks, &context, 1);
  259. BlockRandomAccessDenseMatrix expected_lhs(blocks, &context, 1);
  260. Vector actual_rhs(kFBlockSize);
  261. Vector expected_rhs(kFBlockSize);
  262. Vector f_sol(kFBlockSize);
  263. f_sol.setRandom();
  264. Vector actual_e_sol(num_e_blocks * kEBlockSize);
  265. actual_e_sol.setZero();
  266. Vector expected_e_sol(num_e_blocks * kEBlockSize);
  267. expected_e_sol.setZero();
  268. {
  269. LinearSolver::Options linear_solver_options;
  270. linear_solver_options.e_block_size = kEBlockSize;
  271. linear_solver_options.row_block_size = kRowBlockSize;
  272. linear_solver_options.f_block_size = kFBlockSize;
  273. linear_solver_options.context = &context;
  274. std::unique_ptr<SchurEliminatorBase> eliminator(
  275. SchurEliminatorBase::Create(linear_solver_options));
  276. eliminator->Init(num_e_blocks, true, matrix.block_structure());
  277. eliminator->Eliminate(BlockSparseMatrixData(matrix),
  278. b.data(),
  279. diagonal.data(),
  280. &expected_lhs,
  281. expected_rhs.data());
  282. eliminator->BackSubstitute(BlockSparseMatrixData(matrix),
  283. b.data(),
  284. diagonal.data(),
  285. f_sol.data(),
  286. actual_e_sol.data());
  287. }
  288. {
  289. SchurEliminatorForOneFBlock<2, 3, 6> eliminator;
  290. eliminator.Init(num_e_blocks, true, matrix.block_structure());
  291. eliminator.Eliminate(BlockSparseMatrixData(matrix),
  292. b.data(),
  293. diagonal.data(),
  294. &actual_lhs,
  295. actual_rhs.data());
  296. eliminator.BackSubstitute(BlockSparseMatrixData(matrix),
  297. b.data(),
  298. diagonal.data(),
  299. f_sol.data(),
  300. expected_e_sol.data());
  301. }
  302. ConstMatrixRef actual_lhsref(
  303. actual_lhs.values(), actual_lhs.num_cols(), actual_lhs.num_cols());
  304. ConstMatrixRef expected_lhsref(
  305. expected_lhs.values(), actual_lhs.num_cols(), actual_lhs.num_cols());
  306. EXPECT_NEAR((actual_lhsref - expected_lhsref).norm() / expected_lhsref.norm(),
  307. 0.0,
  308. 1e-12)
  309. << "expected: \n"
  310. << expected_lhsref << "\nactual: \n"
  311. << actual_lhsref;
  312. EXPECT_NEAR(
  313. (actual_rhs - expected_rhs).norm() / expected_rhs.norm(), 0.0, 1e-12)
  314. << "expected: \n"
  315. << expected_rhs << "\nactual: \n"
  316. << actual_rhs;
  317. EXPECT_NEAR((actual_e_sol - expected_e_sol).norm() / expected_e_sol.norm(),
  318. 0.0,
  319. 1e-12)
  320. << "expected: \n"
  321. << expected_e_sol << "\nactual: \n"
  322. << actual_e_sol;
  323. }
  324. } // namespace ceres::internal