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