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