blob: 390f0818140973c92308a2985d2eb33ff7d19cb2 [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#ifndef CERES_NO_SUITESPARSE
32
33#include "ceres/visibility_based_preconditioner.h"
34
35#include <glog/logging.h>
36#include "ceres/file.h"
37#include "gtest/gtest.h"
38#include "Eigen/Dense"
39#include "ceres/block_random_access_dense_matrix.h"
40#include "ceres/block_random_access_sparse_matrix.h"
41#include "ceres/block_sparse_matrix.h"
42#include "ceres/casts.h"
43#include "ceres/collections_port.h"
44#include "ceres/linear_least_squares_problems.h"
45#include "ceres/schur_eliminator.h"
46#include "ceres/stringprintf.h"
47#include "ceres/internal/eigen.h"
48#include "ceres/internal/scoped_ptr.h"
49#include "ceres/types.h"
50
51DECLARE_string(test_srcdir);
52
53namespace ceres {
54namespace internal {
55
56using testing::AssertionResult;
57using testing::AssertionSuccess;
58using testing::AssertionFailure;
59
60static const double kTolerance = 1e-12;
61
62class VisibilityBasedPreconditionerTest : public ::testing::Test {
63 public:
64 static const int kCameraSize = 9;
65
66 protected:
67 void SetUp() {
68 string input_file =
69 JoinPath(FLAGS_test_srcdir,
70 "problem-6-1384-000.lsqp"); // NOLINT
71
72 scoped_ptr<LinearLeastSquaresProblem> problem(
73 CHECK_NOTNULL(CreateLinearLeastSquaresProblemFromFile(input_file)));
74 A_.reset(down_cast<BlockSparseMatrix*>(problem->A.release()));
75 b_.reset(problem->b.release());
76 D_.reset(problem->D.release());
77
78 const CompressedRowBlockStructure* bs =
79 CHECK_NOTNULL(A_->block_structure());
80 const int num_col_blocks = bs->cols.size();
81
82 num_cols_ = A_->num_cols();
83 num_rows_ = A_->num_rows();
84 num_eliminate_blocks_ = problem->num_eliminate_blocks;
85 num_camera_blocks_ = num_col_blocks - num_eliminate_blocks_;
86 options_.num_eliminate_blocks = num_eliminate_blocks_;
87
88 vector<int> blocks(num_col_blocks - num_eliminate_blocks_, 0);
89 for (int i = num_eliminate_blocks_; i < num_col_blocks; ++i) {
90 blocks[i - num_eliminate_blocks_] = bs->cols[i].size;
91 }
92
93 // The input matrix is a real jacobian and fairly poorly
94 // conditioned. Setting D to a large constant makes the normal
95 // equations better conditioned and makes the tests below better
96 // conditioned.
97 VectorRef(D_.get(), num_cols_).setConstant(10.0);
98
99 schur_complement_.reset(new BlockRandomAccessDenseMatrix(blocks));
100 Vector rhs(schur_complement_->num_rows());
101
102 scoped_ptr<SchurEliminatorBase> eliminator;
103 eliminator.reset(SchurEliminatorBase::Create(options_));
104 eliminator->Init(num_eliminate_blocks_, bs);
105 eliminator->Eliminate(A_.get(), b_.get(), D_.get(),
106 schur_complement_.get(), rhs.data());
107 }
108
109
110 AssertionResult IsSparsityStructureValid() {
111 preconditioner_->InitStorage(*A_->block_structure());
112 const HashSet<pair<int, int> >& cluster_pairs = get_cluster_pairs();
113 const vector<int>& cluster_membership = get_cluster_membership();
114
115 for (int i = 0; i < num_camera_blocks_; ++i) {
116 for (int j = i; j < num_camera_blocks_; ++j) {
117 if (cluster_pairs.count(make_pair(cluster_membership[i],
118 cluster_membership[j]))) {
119 if (!IsBlockPairInPreconditioner(i, j)) {
120 return AssertionFailure()
121 << "block pair (" << i << "," << j << "missing";
122 }
123 } else {
124 if (IsBlockPairInPreconditioner(i, j)) {
125 return AssertionFailure()
126 << "block pair (" << i << "," << j << "should not be present";
127 }
128 }
129 }
130 }
131 return AssertionSuccess();
132 }
133
134 AssertionResult PreconditionerValuesMatch() {
135 preconditioner_->Compute(*A_, D_.get());
136 const HashSet<pair<int, int> >& cluster_pairs = get_cluster_pairs();
137 const BlockRandomAccessSparseMatrix* m = get_m();
138 Matrix preconditioner_matrix;
139 m->matrix()->ToDenseMatrix(&preconditioner_matrix);
140 ConstMatrixRef full_schur_complement(schur_complement_->values(),
141 m->num_rows(),
142 m->num_rows());
143 const int num_clusters = get_num_clusters();
144 const int kDiagonalBlockSize =
145 kCameraSize * num_camera_blocks_ / num_clusters;
146
147 for (int i = 0; i < num_clusters; ++i) {
148 for (int j = i; j < num_clusters; ++j) {
149 double diff = 0.0;
150 if (cluster_pairs.count(make_pair(i, j))) {
151 diff =
152 (preconditioner_matrix.block(kDiagonalBlockSize * i,
153 kDiagonalBlockSize * j,
154 kDiagonalBlockSize,
155 kDiagonalBlockSize) -
156 full_schur_complement.block(kDiagonalBlockSize * i,
157 kDiagonalBlockSize * j,
158 kDiagonalBlockSize,
159 kDiagonalBlockSize)).norm();
160 } else {
161 diff = preconditioner_matrix.block(kDiagonalBlockSize * i,
162 kDiagonalBlockSize * j,
163 kDiagonalBlockSize,
164 kDiagonalBlockSize).norm();
165 }
166 if (diff > kTolerance) {
167 return AssertionFailure()
168 << "Preconditioner block " << i << " " << j << " differs "
169 << "from expected value by " << diff;
170 }
171 }
172 }
173 return AssertionSuccess();
174 }
175
176 // Accessors
177 int get_num_blocks() { return preconditioner_->num_blocks_; }
178
179 int get_num_clusters() { return preconditioner_->num_clusters_; }
180 int* get_mutable_num_clusters() { return &preconditioner_->num_clusters_; }
181
182 const vector<int>& get_block_size() {
183 return preconditioner_->block_size_; }
184
185 vector<int>* get_mutable_block_size() {
186 return &preconditioner_->block_size_; }
187
188 const vector<int>& get_cluster_membership() {
189 return preconditioner_->cluster_membership_;
190 }
191
192 vector<int>* get_mutable_cluster_membership() {
193 return &preconditioner_->cluster_membership_;
194 }
195
196 const set<pair<int, int> >& get_block_pairs() {
197 return preconditioner_->block_pairs_;
198 }
199
200 set<pair<int, int> >* get_mutable_block_pairs() {
201 return &preconditioner_->block_pairs_;
202 }
203
204 const HashSet<pair<int, int> >& get_cluster_pairs() {
205 return preconditioner_->cluster_pairs_;
206 }
207
208 HashSet<pair<int, int> >* get_mutable_cluster_pairs() {
209 return &preconditioner_->cluster_pairs_;
210 }
211
212 bool IsBlockPairInPreconditioner(const int block1, const int block2) {
213 return preconditioner_->IsBlockPairInPreconditioner(block1, block2);
214 }
215
216 bool IsBlockPairOffDiagonal(const int block1, const int block2) {
217 return preconditioner_->IsBlockPairOffDiagonal(block1, block2);
218 }
219
220 const BlockRandomAccessSparseMatrix* get_m() {
221 return preconditioner_->m_.get();
222 }
223
224 int num_rows_;
225 int num_cols_;
226 int num_eliminate_blocks_;
227 int num_camera_blocks_;
228
229 scoped_ptr<BlockSparseMatrix> A_;
230 scoped_array<double> b_;
231 scoped_array<double> D_;
232
233 LinearSolver::Options options_;
234 scoped_ptr<VisibilityBasedPreconditioner> preconditioner_;
235 scoped_ptr<BlockRandomAccessDenseMatrix> schur_complement_;
236};
237
238#ifndef CERES_DONT_HAVE_PROTOCOL_BUFFERS
239TEST_F(VisibilityBasedPreconditionerTest, SchurJacobiStructure) {
240 options_.preconditioner_type = SCHUR_JACOBI;
241 preconditioner_.reset(
242 new VisibilityBasedPreconditioner(*A_->block_structure(), options_));
243 EXPECT_EQ(get_num_blocks(), num_camera_blocks_);
244 EXPECT_EQ(get_num_clusters(), num_camera_blocks_);
245 for (int i = 0; i < num_camera_blocks_; ++i) {
246 for (int j = 0; j < num_camera_blocks_; ++j) {
247 const string msg = StringPrintf("Camera pair: %d %d", i, j);
248 SCOPED_TRACE(msg);
249 if (i == j) {
250 EXPECT_TRUE(IsBlockPairInPreconditioner(i, j));
251 EXPECT_FALSE(IsBlockPairOffDiagonal(i, j));
252 } else {
253 EXPECT_FALSE(IsBlockPairInPreconditioner(i, j));
254 EXPECT_TRUE(IsBlockPairOffDiagonal(i, j));
255 }
256 }
257 }
258}
259
260TEST_F(VisibilityBasedPreconditionerTest, OneClusterClusterJacobi) {
261 options_.preconditioner_type = CLUSTER_JACOBI;
262 preconditioner_.reset(
263 new VisibilityBasedPreconditioner(*A_->block_structure(), options_));
264
265 // Override the clustering to be a single clustering containing all
266 // the cameras.
267 vector<int>& cluster_membership = *get_mutable_cluster_membership();
268 for (int i = 0; i < num_camera_blocks_; ++i) {
269 cluster_membership[i] = 0;
270 }
271
272 *get_mutable_num_clusters() = 1;
273
274 HashSet<pair<int, int> >& cluster_pairs = *get_mutable_cluster_pairs();
275 cluster_pairs.clear();
276 cluster_pairs.insert(make_pair(0, 0));
277
278 EXPECT_TRUE(IsSparsityStructureValid());
279 EXPECT_TRUE(PreconditionerValuesMatch());
280
281 // Multiplication by the inverse of the preconditioner.
282 const int num_rows = schur_complement_->num_rows();
283 ConstMatrixRef full_schur_complement(schur_complement_->values(),
284 num_rows,
285 num_rows);
286 Vector x(num_rows);
287 Vector y(num_rows);
288 Vector z(num_rows);
289
290 for (int i = 0; i < num_rows; ++i) {
291 x.setZero();
292 y.setZero();
293 z.setZero();
294 x[i] = 1.0;
295 preconditioner_->RightMultiply(x.data(), y.data());
296 z = full_schur_complement
297 .selfadjointView<Eigen::Upper>()
298 .ldlt().solve(x);
299 double max_relative_difference =
300 ((y - z).array() / z.array()).matrix().lpNorm<Eigen::Infinity>();
301 EXPECT_NEAR(max_relative_difference, 0.0, kTolerance);
302 }
303}
304
305
306
307TEST_F(VisibilityBasedPreconditionerTest, ClusterJacobi) {
308 options_.preconditioner_type = CLUSTER_JACOBI;
309 preconditioner_.reset(
310 new VisibilityBasedPreconditioner(*A_->block_structure(), options_));
311
312 // Override the clustering to be equal number of cameras.
313 vector<int>& cluster_membership = *get_mutable_cluster_membership();
314 cluster_membership.resize(num_camera_blocks_);
315 static const int kNumClusters = 3;
316
317 for (int i = 0; i < num_camera_blocks_; ++i) {
318 cluster_membership[i] = (i * kNumClusters) / num_camera_blocks_;
319 }
320 *get_mutable_num_clusters() = kNumClusters;
321
322 HashSet<pair<int, int> >& cluster_pairs = *get_mutable_cluster_pairs();
323 cluster_pairs.clear();
324 for (int i = 0; i < kNumClusters; ++i) {
325 cluster_pairs.insert(make_pair(i, i));
326 }
327
328 EXPECT_TRUE(IsSparsityStructureValid());
329 EXPECT_TRUE(PreconditionerValuesMatch());
330}
331
332
333TEST_F(VisibilityBasedPreconditionerTest, ClusterTridiagonal) {
334 options_.preconditioner_type = CLUSTER_TRIDIAGONAL;
335 preconditioner_.reset(
336 new VisibilityBasedPreconditioner(*A_->block_structure(), options_));
337 static const int kNumClusters = 3;
338
339 // Override the clustering to be 3 clusters.
340 vector<int>& cluster_membership = *get_mutable_cluster_membership();
341 cluster_membership.resize(num_camera_blocks_);
342 for (int i = 0; i < num_camera_blocks_; ++i) {
343 cluster_membership[i] = (i * kNumClusters) / num_camera_blocks_;
344 }
345 *get_mutable_num_clusters() = kNumClusters;
346
347 // Spanning forest has structure 0-1 2
348 HashSet<pair<int, int> >& cluster_pairs = *get_mutable_cluster_pairs();
349 cluster_pairs.clear();
350 for (int i = 0; i < kNumClusters; ++i) {
351 cluster_pairs.insert(make_pair(i, i));
352 }
353 cluster_pairs.insert(make_pair(0, 1));
354
355 EXPECT_TRUE(IsSparsityStructureValid());
356 EXPECT_TRUE(PreconditionerValuesMatch());
357}
358#endif // CERES_DONT_HAVE_PROTOCOL_BUFFERS
359
360} // namespace internal
361} // namespace ceres
362
363#endif // CERES_NO_SUITESPARSE