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