blob: 0285680d722fbbafbd50f2fa901042270604ce45 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2017 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 <memory>
32#include "ceres/subset_preconditioner.h"
33#include "Eigen/Dense"
34#include "Eigen/SparseCore"
35#include "ceres/block_sparse_matrix.h"
36#include "ceres/compressed_row_sparse_matrix.h"
37#include "ceres/inner_product_computer.h"
38#include "ceres/internal/eigen.h"
39#include "glog/logging.h"
40#include "gtest/gtest.h"
41
42namespace ceres {
43namespace internal {
44
45// TODO(sameeragarwal): Refactor the following two functions out of
46// here and sparse_cholesky_test.cc into a more suitable place.
47template <int UpLoType>
48bool SolveLinearSystemUsingEigen(const Matrix& lhs,
49 const Vector rhs,
50 Vector* solution) {
51 Eigen::LLT<Matrix, UpLoType> llt = lhs.selfadjointView<UpLoType>().llt();
52 if (llt.info() != Eigen::Success) {
53 return false;
54 }
55 *solution = llt.solve(rhs);
56 return (llt.info() == Eigen::Success);
57}
58
59// Use Eigen's Dense Cholesky solver to compute the solution to a
60// sparse linear system.
61bool ComputeExpectedSolution(const CompressedRowSparseMatrix& lhs,
62 const Vector& rhs,
63 Vector* solution) {
64 Matrix dense_triangular_lhs;
65 lhs.ToDenseMatrix(&dense_triangular_lhs);
66 if (lhs.storage_type() == CompressedRowSparseMatrix::UPPER_TRIANGULAR) {
67 Matrix full_lhs = dense_triangular_lhs.selfadjointView<Eigen::Upper>();
68 return SolveLinearSystemUsingEigen<Eigen::Upper>(full_lhs, rhs, solution);
69 }
70 return SolveLinearSystemUsingEigen<Eigen::Lower>(
71 dense_triangular_lhs, rhs, solution);
72}
73
74typedef ::testing::tuple<SparseLinearAlgebraLibraryType, bool> Param;
75
76std::string ParamInfoToString(testing::TestParamInfo<Param> info) {
77 Param param = info.param;
78 std::stringstream ss;
79 ss << SparseLinearAlgebraLibraryTypeToString(::testing::get<0>(param)) << "_"
80 << (::testing::get<1>(param) ? "Diagonal" : "NoDiagonal");
81 return ss.str();
82}
83
84class SubsetPreconditionerTest : public ::testing::TestWithParam<Param> {
85 protected:
86 virtual void SetUp() {
87 BlockSparseMatrix::RandomMatrixOptions options;
88 options.num_col_blocks = 4;
89 options.min_col_block_size = 1;
90 options.max_col_block_size = 4;
91 options.num_row_blocks = 8;
92 options.min_row_block_size = 1;
93 options.max_row_block_size = 4;
94 options.block_density = 0.9;
95
96 m_.reset(BlockSparseMatrix::CreateRandomMatrix(options));
97 start_row_block_ = m_->block_structure()->rows.size();
98
99 // Ensure that the bottom part of the matrix has the same column
100 // block structure.
101 options.col_blocks = m_->block_structure()->cols;
102 b_.reset(BlockSparseMatrix::CreateRandomMatrix(options));
103 m_->AppendRows(*b_);
104
105 // Create a Identity block diagonal matrix with the same column
106 // block structure.
107 diagonal_ = Vector::Ones(m_->num_cols());
108 block_diagonal_.reset(BlockSparseMatrix::CreateDiagonalMatrix(
109 diagonal_.data(), b_->block_structure()->cols));
110
111 // Unconditionally add the block diagonal to the matrix b_,
112 // because either it is either part of b_ to make it full rank, or
113 // we pass the same diagonal matrix later as the parameter D. In
114 // either case the preconditioner matrix is b_' b + D'D.
115 b_->AppendRows(*block_diagonal_);
116 inner_product_computer_.reset(InnerProductComputer::Create(
117 *b_, CompressedRowSparseMatrix::UPPER_TRIANGULAR));
118 inner_product_computer_->Compute();
119 }
120
121 std::unique_ptr<BlockSparseMatrix> m_;
122 std::unique_ptr<BlockSparseMatrix> b_;
123 std::unique_ptr<BlockSparseMatrix> block_diagonal_;
124 std::unique_ptr<InnerProductComputer> inner_product_computer_;
125 std::unique_ptr<Preconditioner> preconditioner_;
126 Vector diagonal_;
127 int start_row_block_;
128};
129
130TEST_P(SubsetPreconditionerTest, foo) {
131 Param param = GetParam();
132 Preconditioner::Options options;
133 options.subset_preconditioner_start_row_block = start_row_block_;
134 options.sparse_linear_algebra_library_type = ::testing::get<0>(param);
135 preconditioner_.reset(new SubsetPreconditioner(options, *m_));
136
137 const bool with_diagonal = ::testing::get<1>(param);
138 if (!with_diagonal) {
139 m_->AppendRows(*block_diagonal_);
140 }
141
142 EXPECT_TRUE(
143 preconditioner_->Update(*m_, with_diagonal ? diagonal_.data() : NULL));
144
145 // Repeatedly apply the preconditioner to random vectors and check
146 // that the preconditioned value is the same as one obtained by
147 // solving the linear system directly.
148 for (int i = 0; i < 5; ++i) {
149 CompressedRowSparseMatrix* lhs = inner_product_computer_->mutable_result();
150 Vector rhs = Vector::Random(lhs->num_rows());
151 Vector expected(lhs->num_rows());
152 EXPECT_TRUE(ComputeExpectedSolution(*lhs, rhs, &expected));
153
154 Vector actual(lhs->num_rows());
155 preconditioner_->RightMultiply(rhs.data(), actual.data());
156
157 Matrix eigen_lhs;
158 lhs->ToDenseMatrix(&eigen_lhs);
159 EXPECT_NEAR((actual - expected).norm() / actual.norm(),
160 0.0,
161 std::numeric_limits<double>::epsilon() * 10)
162 << "\n"
163 << eigen_lhs << "\n"
164 << expected.transpose() << "\n"
165 << actual.transpose();
166 }
167}
168
169#ifndef CERES_NO_SUITESPARSE
170INSTANTIATE_TEST_CASE_P(SubsetPreconditionerWithSuiteSparse,
171 SubsetPreconditionerTest,
172 ::testing::Combine(::testing::Values(SUITE_SPARSE),
173 ::testing::Values(true, false)),
174 ParamInfoToString);
175#endif
176
177#ifndef CERES_NO_CXSPARSE
178INSTANTIATE_TEST_CASE_P(SubsetPreconditionerWithCXSparse,
179 SubsetPreconditionerTest,
180 ::testing::Combine(::testing::Values(CX_SPARSE),
181 ::testing::Values(true, false)),
182 ParamInfoToString);
183#endif
184
185#ifndef CERES_NO_ACCELERATE_SPARSE
186INSTANTIATE_TEST_CASE_P(SubsetPreconditionerWithAccelerateSparse,
187 SubsetPreconditionerTest,
188 ::testing::Combine(::testing::Values(ACCELERATE_SPARSE),
189 ::testing::Values(true, false)),
190 ParamInfoToString);
191#endif
192
193#ifdef CERES_USE_EIGEN_SPARSE
194INSTANTIATE_TEST_CASE_P(SubsetPreconditionerWithEigenSparse,
195 SubsetPreconditionerTest,
196 ::testing::Combine(::testing::Values(EIGEN_SPARSE),
197 ::testing::Values(true, false)),
198 ParamInfoToString);
199#endif
200
201} // namespace internal
202} // namespace ceres