blob: 0ac435475f309c000b4bdd2d57b8622a9a74a513 [file] [log] [blame]
Keir Mierle8ebb0732012-04-30 23:09:08 -07001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2010, 2011, 2012 Google Inc. All rights reserved.
3// http://code.google.com/p/ceres-solver/
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 <glog/logging.h>
34#include "ceres/file.h"
35#include "gtest/gtest.h"
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/detect_structure.h"
41#include "ceres/linear_least_squares_problems.h"
42#include "ceres/triplet_sparse_matrix.h"
43#include "ceres/internal/eigen.h"
44#include "ceres/internal/scoped_ptr.h"
45#include "ceres/types.h"
46
47// TODO(sameeragarwal): Reduce the size of these tests and redo the
48// parameterization to be more efficient.
49
50DECLARE_string(test_srcdir);
51
52namespace ceres {
53namespace internal {
54
55class SchurEliminatorTest : public ::testing::Test {
56 protected:
57 void SetUpFromId(int id) {
58 scoped_ptr<LinearLeastSquaresProblem>
59 problem(CreateLinearLeastSquaresProblemFromId(id));
60 CHECK_NOTNULL(problem.get());
61 SetupHelper(problem.get());
62 }
63
64 void SetUpFromFilename(const string& filename) {
65 scoped_ptr<LinearLeastSquaresProblem>
66 problem(CreateLinearLeastSquaresProblemFromFile(filename));
67 CHECK_NOTNULL(problem.get());
68 SetupHelper(problem.get());
69 }
70
71 void SetupHelper(LinearLeastSquaresProblem* problem) {
72 A.reset(down_cast<BlockSparseMatrix*>(problem->A.release()));
73 b.reset(problem->b.release());
74 D.reset(problem->D.release());
75
76 num_eliminate_blocks = problem->num_eliminate_blocks;
77 num_eliminate_cols = 0;
78 const CompressedRowBlockStructure* bs = A->block_structure();
79
80 for (int i = 0; i < num_eliminate_blocks; ++i) {
81 num_eliminate_cols += bs->cols[i].size;
82 }
83 }
84
85 // Compute the golden values for the reduced linear system and the
86 // solution to the linear least squares problem using dense linear
87 // algebra.
88 void ComputeReferenceSolution(const Vector& D) {
89 Matrix J;
90 A->ToDenseMatrix(&J);
91 VectorRef f(b.get(), J.rows());
92
93 Matrix H = (D.cwiseProduct(D)).asDiagonal();
94 H.noalias() += J.transpose() * J;
95
96 const Vector g = J.transpose() * f;
97 const int schur_size = J.cols() - num_eliminate_cols;
98
99 lhs_expected.resize(schur_size, schur_size);
100 lhs_expected.setZero();
101
102 rhs_expected.resize(schur_size);
103 rhs_expected.setZero();
104
105 sol_expected.resize(J.cols());
106 sol_expected.setZero();
107
108 Matrix P = H.block(0, 0, num_eliminate_cols, num_eliminate_cols);
109 Matrix Q = H.block(0,
110 num_eliminate_cols,
111 num_eliminate_cols,
112 schur_size);
113 Matrix R = H.block(num_eliminate_cols,
114 num_eliminate_cols,
115 schur_size,
116 schur_size);
117 int row = 0;
118 const CompressedRowBlockStructure* bs = A->block_structure();
119 for (int i = 0; i < num_eliminate_blocks; ++i) {
120 const int block_size = bs->cols[i].size;
121 P.block(row, row, block_size, block_size) =
122 P
123 .block(row, row, block_size, block_size)
124 .ldlt()
125 .solve(Matrix::Identity(block_size, block_size));
126 row += block_size;
127 }
128
129 lhs_expected
130 .triangularView<Eigen::Upper>() = R - Q.transpose() * P * Q;
131 rhs_expected =
132 g.tail(schur_size) - Q.transpose() * P * g.head(num_eliminate_cols);
133 sol_expected = H.ldlt().solve(g);
134 }
135
136 void EliminateSolveAndCompare(const VectorRef& diagonal,
137 bool use_static_structure,
138 const double relative_tolerance) {
139 const CompressedRowBlockStructure* bs = A->block_structure();
140 const int num_col_blocks = bs->cols.size();
141 vector<int> blocks(num_col_blocks - num_eliminate_blocks, 0);
142 for (int i = num_eliminate_blocks; i < num_col_blocks; ++i) {
143 blocks[i - num_eliminate_blocks] = bs->cols[i].size;
144 }
145
146 BlockRandomAccessDenseMatrix lhs(blocks);
147
148 const int num_cols = A->num_cols();
149 const int schur_size = lhs.num_rows();
150
151 Vector rhs(schur_size);
152
153 LinearSolver::Options options;
154 options.num_eliminate_blocks = num_eliminate_blocks;
155 if (use_static_structure) {
156 DetectStructure(*bs,
157 options.num_eliminate_blocks,
158 &options.row_block_size,
159 &options.e_block_size,
160 &options.f_block_size);
161 }
162
163 scoped_ptr<SchurEliminatorBase> eliminator;
164 eliminator.reset(SchurEliminatorBase::Create(options));
165 eliminator->Init(num_eliminate_blocks, A->block_structure());
166 eliminator->Eliminate(A.get(), b.get(), diagonal.data(), &lhs, rhs.data());
167
168 MatrixRef lhs_ref(lhs.mutable_values(), lhs.num_rows(), lhs.num_cols());
169 Vector reduced_sol =
170 lhs_ref
171 .selfadjointView<Eigen::Upper>()
172 .ldlt()
173 .solve(rhs);
174
175 // Solution to the linear least squares problem.
176 Vector sol(num_cols);
177 sol.setZero();
178 sol.tail(schur_size) = reduced_sol;
179 eliminator->BackSubstitute(A.get(),
180 b.get(),
181 diagonal.data(),
182 reduced_sol.data(),
183 sol.data());
184
185 Matrix delta = (lhs_ref - lhs_expected).selfadjointView<Eigen::Upper>();
186 double diff = delta.norm();
187 EXPECT_NEAR(diff / lhs_expected.norm(), 0.0, relative_tolerance);
188 EXPECT_NEAR((rhs - rhs_expected).norm() / rhs_expected.norm(), 0.0,
189 relative_tolerance);
190 EXPECT_NEAR((sol - sol_expected).norm() / sol_expected.norm(), 0.0,
191 relative_tolerance);
192 }
193
194 scoped_ptr<BlockSparseMatrix> A;
195 scoped_array<double> b;
196 scoped_array<double> D;
197 int num_eliminate_blocks;
198 int num_eliminate_cols;
199
200 Matrix lhs_expected;
201 Vector rhs_expected;
202 Vector sol_expected;
203};
204
205TEST_F(SchurEliminatorTest, ScalarProblem) {
206 SetUpFromId(2);
207 Vector zero(A->num_cols());
208 zero.setZero();
209
210 ComputeReferenceSolution(VectorRef(zero.data(), A->num_cols()));
211 EliminateSolveAndCompare(VectorRef(zero.data(), A->num_cols()), true, 1e-14);
212 EliminateSolveAndCompare(VectorRef(zero.data(), A->num_cols()), false, 1e-14);
213
214 ComputeReferenceSolution(VectorRef(D.get(), A->num_cols()));
215 EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), true, 1e-14);
216 EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), false, 1e-14);
217}
218
219#ifndef CERES_DONT_HAVE_PROTOCOL_BUFFERS
220TEST_F(SchurEliminatorTest, BlockProblem) {
221 const string input_file =
222 JoinPath(FLAGS_test_srcdir,
223 "problem-6-1384-000.lsqp"); // NOLINT
224
225 SetUpFromFilename(input_file);
226 ComputeReferenceSolution(VectorRef(D.get(), A->num_cols()));
227 EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), true, 1e-10);
228 EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), false, 1e-10);
229}
230#endif // CERES_DONT_HAVE_PROTOCOL_BUFFERS
231
232} // namespace internal
233} // namespace ceres