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