blob: 4c1f6d8f4905ee5d3b34a0c9027d789588038bb7 [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 "ceres/sparse_cholesky.h"
32
33#include <memory>
34#include <numeric>
35#include <vector>
36
37#include "Eigen/Dense"
38#include "Eigen/SparseCore"
39#include "ceres/block_sparse_matrix.h"
40#include "ceres/compressed_row_sparse_matrix.h"
41#include "ceres/inner_product_computer.h"
42#include "ceres/internal/eigen.h"
43#include "ceres/iterative_refiner.h"
44#include "ceres/random.h"
45#include "glog/logging.h"
46#include "gmock/gmock.h"
47#include "gtest/gtest.h"
48
49namespace ceres {
50namespace internal {
51
52BlockSparseMatrix* CreateRandomFullRankMatrix(const int num_col_blocks,
53 const int min_col_block_size,
54 const int max_col_block_size,
55 const double block_density) {
56 // Create a random matrix
57 BlockSparseMatrix::RandomMatrixOptions options;
58 options.num_col_blocks = num_col_blocks;
59 options.min_col_block_size = min_col_block_size;
60 options.max_col_block_size = max_col_block_size;
61
62 options.num_row_blocks = 2 * num_col_blocks;
63 options.min_row_block_size = 1;
64 options.max_row_block_size = max_col_block_size;
65 options.block_density = block_density;
66 std::unique_ptr<BlockSparseMatrix> random_matrix(
67 BlockSparseMatrix::CreateRandomMatrix(options));
68
69 // Add a diagonal block sparse matrix to make it full rank.
70 Vector diagonal = Vector::Ones(random_matrix->num_cols());
71 std::unique_ptr<BlockSparseMatrix> block_diagonal(
72 BlockSparseMatrix::CreateDiagonalMatrix(
73 diagonal.data(), random_matrix->block_structure()->cols));
74 random_matrix->AppendRows(*block_diagonal);
75 return random_matrix.release();
76}
77
78bool ComputeExpectedSolution(const CompressedRowSparseMatrix& lhs,
79 const Vector& rhs,
80 Vector* solution) {
81 Matrix eigen_lhs;
82 lhs.ToDenseMatrix(&eigen_lhs);
83 if (lhs.storage_type() == CompressedRowSparseMatrix::UPPER_TRIANGULAR) {
84 Matrix full_lhs = eigen_lhs.selfadjointView<Eigen::Upper>();
85 Eigen::LLT<Matrix, Eigen::Upper> llt =
86 eigen_lhs.selfadjointView<Eigen::Upper>().llt();
87 if (llt.info() != Eigen::Success) {
88 return false;
89 }
90 *solution = llt.solve(rhs);
91 return (llt.info() == Eigen::Success);
92 }
93
94 Matrix full_lhs = eigen_lhs.selfadjointView<Eigen::Lower>();
95 Eigen::LLT<Matrix, Eigen::Lower> llt =
96 eigen_lhs.selfadjointView<Eigen::Lower>().llt();
97 if (llt.info() != Eigen::Success) {
98 return false;
99 }
100 *solution = llt.solve(rhs);
101 return (llt.info() == Eigen::Success);
102}
103
104void SparseCholeskySolverUnitTest(
105 const SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
106 const OrderingType ordering_type,
107 const bool use_block_structure,
108 const int num_blocks,
109 const int min_block_size,
110 const int max_block_size,
111 const double block_density) {
112 LinearSolver::Options sparse_cholesky_options;
113 sparse_cholesky_options.sparse_linear_algebra_library_type =
114 sparse_linear_algebra_library_type;
115 sparse_cholesky_options.use_postordering = (ordering_type == AMD);
116 std::unique_ptr<SparseCholesky> sparse_cholesky = SparseCholesky::Create(
117 sparse_cholesky_options);
118 const CompressedRowSparseMatrix::StorageType storage_type =
119 sparse_cholesky->StorageType();
120
121 std::unique_ptr<BlockSparseMatrix> m(CreateRandomFullRankMatrix(
122 num_blocks, min_block_size, max_block_size, block_density));
123 std::unique_ptr<InnerProductComputer> inner_product_computer(
124 InnerProductComputer::Create(*m, storage_type));
125 inner_product_computer->Compute();
126 CompressedRowSparseMatrix* lhs = inner_product_computer->mutable_result();
127
128 if (!use_block_structure) {
129 lhs->mutable_row_blocks()->clear();
130 lhs->mutable_col_blocks()->clear();
131 }
132
133 Vector rhs = Vector::Random(lhs->num_rows());
134 Vector expected(lhs->num_rows());
135 Vector actual(lhs->num_rows());
136
137 EXPECT_TRUE(ComputeExpectedSolution(*lhs, rhs, &expected));
138 std::string message;
139 EXPECT_EQ(sparse_cholesky->FactorAndSolve(
140 lhs, rhs.data(), actual.data(), &message),
141 LINEAR_SOLVER_SUCCESS);
142 Matrix eigen_lhs;
143 lhs->ToDenseMatrix(&eigen_lhs);
144 EXPECT_NEAR((actual - expected).norm() / actual.norm(),
145 0.0,
146 std::numeric_limits<double>::epsilon() * 20)
147 << "\n"
148 << eigen_lhs;
149}
150
151typedef ::testing::tuple<SparseLinearAlgebraLibraryType, OrderingType, bool>
152 Param;
153
154std::string ParamInfoToString(testing::TestParamInfo<Param> info) {
155 Param param = info.param;
156 std::stringstream ss;
157 ss << SparseLinearAlgebraLibraryTypeToString(::testing::get<0>(param)) << "_"
158 << (::testing::get<1>(param) == AMD ? "AMD" : "NATURAL") << "_"
159 << (::testing::get<2>(param) ? "UseBlockStructure" : "NoBlockStructure");
160 return ss.str();
161}
162
163class SparseCholeskyTest : public ::testing::TestWithParam<Param> {};
164
165TEST_P(SparseCholeskyTest, FactorAndSolve) {
166 SetRandomState(2982);
167 const int kMinNumBlocks = 1;
168 const int kMaxNumBlocks = 10;
169 const int kNumTrials = 10;
170 const int kMinBlockSize = 1;
171 const int kMaxBlockSize = 5;
172
173 for (int num_blocks = kMinNumBlocks; num_blocks < kMaxNumBlocks;
174 ++num_blocks) {
175 for (int trial = 0; trial < kNumTrials; ++trial) {
176 const double block_density = std::max(0.1, RandDouble());
177 Param param = GetParam();
178 SparseCholeskySolverUnitTest(::testing::get<0>(param),
179 ::testing::get<1>(param),
180 ::testing::get<2>(param),
181 num_blocks,
182 kMinBlockSize,
183 kMaxBlockSize,
184 block_density);
185 }
186 }
187}
188
189#ifndef CERES_NO_SUITESPARSE
190INSTANTIATE_TEST_CASE_P(SuiteSparseCholesky,
191 SparseCholeskyTest,
192 ::testing::Combine(::testing::Values(SUITE_SPARSE),
193 ::testing::Values(AMD, NATURAL),
194 ::testing::Values(true, false)),
195 ParamInfoToString);
196#endif
197
198#ifndef CERES_NO_CXSPARSE
199INSTANTIATE_TEST_CASE_P(CXSparseCholesky,
200 SparseCholeskyTest,
201 ::testing::Combine(::testing::Values(CX_SPARSE),
202 ::testing::Values(AMD, NATURAL),
203 ::testing::Values(true, false)),
204 ParamInfoToString);
205#endif
206
207#ifndef CERES_NO_ACCELERATE_SPARSE
208INSTANTIATE_TEST_CASE_P(AccelerateSparseCholesky,
209 SparseCholeskyTest,
210 ::testing::Combine(::testing::Values(ACCELERATE_SPARSE),
211 ::testing::Values(AMD, NATURAL),
212 ::testing::Values(true, false)),
213 ParamInfoToString);
214
215INSTANTIATE_TEST_CASE_P(AccelerateSparseCholeskySingle,
216 SparseCholeskyTest,
217 ::testing::Combine(::testing::Values(ACCELERATE_SPARSE),
218 ::testing::Values(AMD, NATURAL),
219 ::testing::Values(true, false)),
220 ParamInfoToString);
221#endif
222
223#ifdef CERES_USE_EIGEN_SPARSE
224INSTANTIATE_TEST_CASE_P(EigenSparseCholesky,
225 SparseCholeskyTest,
226 ::testing::Combine(::testing::Values(EIGEN_SPARSE),
227 ::testing::Values(AMD, NATURAL),
228 ::testing::Values(true, false)),
229 ParamInfoToString);
230
231INSTANTIATE_TEST_CASE_P(EigenSparseCholeskySingle,
232 SparseCholeskyTest,
233 ::testing::Combine(::testing::Values(EIGEN_SPARSE),
234 ::testing::Values(AMD, NATURAL),
235 ::testing::Values(true, false)),
236 ParamInfoToString);
237#endif
238
239class MockSparseCholesky : public SparseCholesky {
240 public:
241 MOCK_CONST_METHOD0(StorageType, CompressedRowSparseMatrix::StorageType());
242 MOCK_METHOD2(Factorize,
243 LinearSolverTerminationType(CompressedRowSparseMatrix* lhs,
244 std::string* message));
245 MOCK_METHOD3(Solve,
246 LinearSolverTerminationType(const double* rhs,
247 double* solution,
248 std::string* message));
249};
250
251class MockIterativeRefiner : public IterativeRefiner {
252 public:
253 MockIterativeRefiner() : IterativeRefiner(1) {}
254 MOCK_METHOD4(Refine,
255 void (const SparseMatrix& lhs,
256 const double* rhs,
257 SparseCholesky* sparse_cholesky,
258 double* solution));
259};
260
261
262using testing::_;
263using testing::Return;
264
265TEST(RefinedSparseCholesky, StorageType) {
266 MockSparseCholesky* mock_sparse_cholesky = new MockSparseCholesky;
267 MockIterativeRefiner* mock_iterative_refiner = new MockIterativeRefiner;
268 EXPECT_CALL(*mock_sparse_cholesky, StorageType())
269 .Times(1)
270 .WillRepeatedly(Return(CompressedRowSparseMatrix::UPPER_TRIANGULAR));
271 EXPECT_CALL(*mock_iterative_refiner, Refine(_, _, _, _))
272 .Times(0);
273 std::unique_ptr<SparseCholesky> sparse_cholesky(mock_sparse_cholesky);
274 std::unique_ptr<IterativeRefiner> iterative_refiner(mock_iterative_refiner);
275 RefinedSparseCholesky refined_sparse_cholesky(std::move(sparse_cholesky),
276 std::move(iterative_refiner));
277 EXPECT_EQ(refined_sparse_cholesky.StorageType(),
278 CompressedRowSparseMatrix::UPPER_TRIANGULAR);
279};
280
281TEST(RefinedSparseCholesky, Factorize) {
282 MockSparseCholesky* mock_sparse_cholesky = new MockSparseCholesky;
283 MockIterativeRefiner* mock_iterative_refiner = new MockIterativeRefiner;
284 EXPECT_CALL(*mock_sparse_cholesky, Factorize(_, _))
285 .Times(1)
286 .WillRepeatedly(Return(LINEAR_SOLVER_SUCCESS));
287 EXPECT_CALL(*mock_iterative_refiner, Refine(_, _, _, _))
288 .Times(0);
289 std::unique_ptr<SparseCholesky> sparse_cholesky(mock_sparse_cholesky);
290 std::unique_ptr<IterativeRefiner> iterative_refiner(mock_iterative_refiner);
291 RefinedSparseCholesky refined_sparse_cholesky(std::move(sparse_cholesky),
292 std::move(iterative_refiner));
293 CompressedRowSparseMatrix m(1, 1, 1);
294 std::string message;
295 EXPECT_EQ(refined_sparse_cholesky.Factorize(&m, &message),
296 LINEAR_SOLVER_SUCCESS);
297};
298
299TEST(RefinedSparseCholesky, FactorAndSolveWithUnsuccessfulFactorization) {
300 MockSparseCholesky* mock_sparse_cholesky = new MockSparseCholesky;
301 MockIterativeRefiner* mock_iterative_refiner = new MockIterativeRefiner;
302 EXPECT_CALL(*mock_sparse_cholesky, Factorize(_, _))
303 .Times(1)
304 .WillRepeatedly(Return(LINEAR_SOLVER_FAILURE));
305 EXPECT_CALL(*mock_sparse_cholesky, Solve(_, _, _))
306 .Times(0);
307 EXPECT_CALL(*mock_iterative_refiner, Refine(_, _, _, _))
308 .Times(0);
309 std::unique_ptr<SparseCholesky> sparse_cholesky(mock_sparse_cholesky);
310 std::unique_ptr<IterativeRefiner> iterative_refiner(mock_iterative_refiner);
311 RefinedSparseCholesky refined_sparse_cholesky(std::move(sparse_cholesky),
312 std::move(iterative_refiner));
313 CompressedRowSparseMatrix m(1, 1, 1);
314 std::string message;
315 double rhs;
316 double solution;
317 EXPECT_EQ(refined_sparse_cholesky.FactorAndSolve(&m, &rhs, &solution, &message),
318 LINEAR_SOLVER_FAILURE);
319};
320
321TEST(RefinedSparseCholesky, FactorAndSolveWithSuccess) {
322 MockSparseCholesky* mock_sparse_cholesky = new MockSparseCholesky;
323 std::unique_ptr<MockIterativeRefiner> mock_iterative_refiner(new MockIterativeRefiner);
324 EXPECT_CALL(*mock_sparse_cholesky, Factorize(_, _))
325 .Times(1)
326 .WillRepeatedly(Return(LINEAR_SOLVER_SUCCESS));
327 EXPECT_CALL(*mock_sparse_cholesky, Solve(_, _, _))
328 .Times(1)
329 .WillRepeatedly(Return(LINEAR_SOLVER_SUCCESS));
330 EXPECT_CALL(*mock_iterative_refiner, Refine(_, _, _, _))
331 .Times(1);
332
333 std::unique_ptr<SparseCholesky> sparse_cholesky(mock_sparse_cholesky);
334 std::unique_ptr<IterativeRefiner> iterative_refiner(std::move(mock_iterative_refiner));
335 RefinedSparseCholesky refined_sparse_cholesky(std::move(sparse_cholesky),
336 std::move(iterative_refiner));
337 CompressedRowSparseMatrix m(1, 1, 1);
338 std::string message;
339 double rhs;
340 double solution;
341 EXPECT_EQ(refined_sparse_cholesky.FactorAndSolve(&m, &rhs, &solution, &message),
342 LINEAR_SOLVER_SUCCESS);
343};
344
345} // namespace internal
346} // namespace ceres