1 // Ceres Solver - A fast non-linear least squares minimizer
2 // Copyright 2015 Google Inc. All rights reserved.
3 // http://ceres-solver.org/
5 // Redistribution and use in source and binary forms, with or without
6 // modification, are permitted provided that the following conditions are met:
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.
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.
29 // Author: sameeragarwal@google.com (Sameer Agarwal)
31 #include "ceres/schur_eliminator.h"
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/detect_structure.h"
38 #include "ceres/internal/eigen.h"
39 #include "ceres/internal/scoped_ptr.h"
40 #include "ceres/linear_least_squares_problems.h"
41 #include "ceres/test_util.h"
42 #include "ceres/triplet_sparse_matrix.h"
43 #include "ceres/types.h"
44 #include "glog/logging.h"
45 #include "gtest/gtest.h"
47 // TODO(sameeragarwal): Reduce the size of these tests and redo the
48 // parameterization to be more efficient.
53 class SchurEliminatorTest : public ::testing::Test {
55 void SetUpFromId(int id) {
56 scoped_ptr<LinearLeastSquaresProblem>
57 problem(CreateLinearLeastSquaresProblemFromId(id));
58 CHECK_NOTNULL(problem.get());
59 SetupHelper(problem.get());
62 void SetupHelper(LinearLeastSquaresProblem* problem) {
63 A.reset(down_cast<BlockSparseMatrix*>(problem->A.release()));
64 b.reset(problem->b.release());
65 D.reset(problem->D.release());
67 num_eliminate_blocks = problem->num_eliminate_blocks;
68 num_eliminate_cols = 0;
69 const CompressedRowBlockStructure* bs = A->block_structure();
71 for (int i = 0; i < num_eliminate_blocks; ++i) {
72 num_eliminate_cols += bs->cols[i].size;
76 // Compute the golden values for the reduced linear system and the
77 // solution to the linear least squares problem using dense linear
79 void ComputeReferenceSolution(const Vector& D) {
82 VectorRef f(b.get(), J.rows());
84 Matrix H = (D.cwiseProduct(D)).asDiagonal();
85 H.noalias() += J.transpose() * J;
87 const Vector g = J.transpose() * f;
88 const int schur_size = J.cols() - num_eliminate_cols;
90 lhs_expected.resize(schur_size, schur_size);
91 lhs_expected.setZero();
93 rhs_expected.resize(schur_size);
94 rhs_expected.setZero();
96 sol_expected.resize(J.cols());
97 sol_expected.setZero();
99 Matrix P = H.block(0, 0, num_eliminate_cols, num_eliminate_cols);
100 Matrix Q = H.block(0,
104 Matrix R = H.block(num_eliminate_cols,
109 const CompressedRowBlockStructure* bs = A->block_structure();
110 for (int i = 0; i < num_eliminate_blocks; ++i) {
111 const int block_size = bs->cols[i].size;
112 P.block(row, row, block_size, block_size) =
114 .block(row, row, block_size, block_size)
116 .solve(Matrix::Identity(block_size, block_size));
121 .triangularView<Eigen::Upper>() = R - Q.transpose() * P * Q;
123 g.tail(schur_size) - Q.transpose() * P * g.head(num_eliminate_cols);
124 sol_expected = H.llt().solve(g);
127 void EliminateSolveAndCompare(const VectorRef& diagonal,
128 bool use_static_structure,
129 const double relative_tolerance) {
130 const CompressedRowBlockStructure* bs = A->block_structure();
131 const int num_col_blocks = bs->cols.size();
132 std::vector<int> blocks(num_col_blocks - num_eliminate_blocks, 0);
133 for (int i = num_eliminate_blocks; i < num_col_blocks; ++i) {
134 blocks[i - num_eliminate_blocks] = bs->cols[i].size;
137 BlockRandomAccessDenseMatrix lhs(blocks);
139 const int num_cols = A->num_cols();
140 const int schur_size = lhs.num_rows();
142 Vector rhs(schur_size);
144 LinearSolver::Options options;
145 options.elimination_groups.push_back(num_eliminate_blocks);
146 if (use_static_structure) {
148 num_eliminate_blocks,
149 &options.row_block_size,
150 &options.e_block_size,
151 &options.f_block_size);
154 scoped_ptr<SchurEliminatorBase> eliminator;
155 eliminator.reset(SchurEliminatorBase::Create(options));
156 const bool kFullRankETE = true;
157 eliminator->Init(num_eliminate_blocks, kFullRankETE, A->block_structure());
158 eliminator->Eliminate(A.get(), b.get(), diagonal.data(), &lhs, rhs.data());
160 MatrixRef lhs_ref(lhs.mutable_values(), lhs.num_rows(), lhs.num_cols());
163 .selfadjointView<Eigen::Upper>()
167 // Solution to the linear least squares problem.
168 Vector sol(num_cols);
170 sol.tail(schur_size) = reduced_sol;
171 eliminator->BackSubstitute(A.get(),
177 Matrix delta = (lhs_ref - lhs_expected).selfadjointView<Eigen::Upper>();
178 double diff = delta.norm();
179 EXPECT_NEAR(diff / lhs_expected.norm(), 0.0, relative_tolerance);
180 EXPECT_NEAR((rhs - rhs_expected).norm() / rhs_expected.norm(), 0.0,
182 EXPECT_NEAR((sol - sol_expected).norm() / sol_expected.norm(), 0.0,
186 scoped_ptr<BlockSparseMatrix> A;
187 scoped_array<double> b;
188 scoped_array<double> D;
189 int num_eliminate_blocks;
190 int num_eliminate_cols;
197 TEST_F(SchurEliminatorTest, ScalarProblemNoRegularization) {
199 Vector zero(A->num_cols());
202 ComputeReferenceSolution(VectorRef(zero.data(), A->num_cols()));
203 EliminateSolveAndCompare(VectorRef(zero.data(), A->num_cols()), true, 1e-14);
204 EliminateSolveAndCompare(VectorRef(zero.data(), A->num_cols()), false, 1e-14);
207 TEST_F(SchurEliminatorTest, ScalarProblemWithRegularization) {
209 ComputeReferenceSolution(VectorRef(D.get(), A->num_cols()));
210 EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), true, 1e-14);
211 EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), false, 1e-14);
214 TEST_F(SchurEliminatorTest, VaryingFBlockSizeWithStaticStructure) {
216 ComputeReferenceSolution(VectorRef(D.get(), A->num_cols()));
217 EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), true, 1e-14);
220 TEST_F(SchurEliminatorTest, VaryingFBlockSizeWithoutStaticStructure) {
222 ComputeReferenceSolution(VectorRef(D.get(), A->num_cols()));
223 EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), false, 1e-14);
226 } // namespace internal