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