1 // Ceres Solver - A fast non-linear least squares minimizer
2 // Copyright 2017 Google Inc. All rights reserved.
3 // http://ceres-solver.org/
5 // Redistribution and use in source and binary forms, with or without
6 // modification, are permitted provided that the following conditions are met:
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.
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.
29 // Author: sameeragarwal@google.com (Sameer Agarwal)
31 #include "ceres/block_sparse_matrix.h"
32 #include "ceres/casts.h"
33 #include "ceres/internal/scoped_ptr.h"
34 #include "ceres/linear_least_squares_problems.h"
35 #include "ceres/linear_solver.h"
36 #include "ceres/triplet_sparse_matrix.h"
37 #include "ceres/types.h"
38 #include "glog/logging.h"
39 #include "gtest/gtest.h"
41 #include "Eigen/Cholesky"
46 // TODO(sameeragarwal): These tests needs to be re-written, since
47 // SparseNormalCholeskySolver is a composition of two classes now,
48 // InnerProductComputer and SparseCholesky.
50 // So the test should exercise the composition, rather than the
51 // numerics of the solver, which are well covered by tests for those
53 class SparseNormalCholeskySolverTest : public ::testing::Test {
55 virtual void SetUp() {
56 scoped_ptr<LinearLeastSquaresProblem> problem(
57 CreateLinearLeastSquaresProblemFromId(2));
59 CHECK_NOTNULL(problem.get());
60 A_.reset(down_cast<BlockSparseMatrix*>(problem->A.release()));
61 b_.reset(problem->b.release());
62 D_.reset(problem->D.release());
65 void TestSolver(const LinearSolver::Options& options, double* D) {
67 A_->ToDenseMatrix(&dense_A);
68 Matrix lhs = dense_A.transpose() * dense_A;
70 lhs += (ConstVectorRef(D, A_->num_cols()).array() *
71 ConstVectorRef(D, A_->num_cols()).array())
76 Vector rhs(A_->num_cols());
78 A_->LeftMultiply(b_.get(), rhs.data());
79 Vector expected_solution = lhs.llt().solve(rhs);
81 scoped_ptr<LinearSolver> solver(LinearSolver::Create(options));
82 LinearSolver::PerSolveOptions per_solve_options;
83 per_solve_options.D = D;
84 Vector actual_solution(A_->num_cols());
85 LinearSolver::Summary summary;
86 summary = solver->Solve(
87 A_.get(), b_.get(), per_solve_options, actual_solution.data());
89 EXPECT_EQ(summary.termination_type, LINEAR_SOLVER_SUCCESS);
91 for (int i = 0; i < A_->num_cols(); ++i) {
92 EXPECT_NEAR(expected_solution(i), actual_solution(i), 1e-8)
93 << "\nExpected: " << expected_solution.transpose()
94 << "\nActual: " << actual_solution.transpose();
98 void TestSolver(const LinearSolver::Options& options) {
99 TestSolver(options, NULL);
100 TestSolver(options, D_.get());
103 scoped_ptr<BlockSparseMatrix> A_;
104 scoped_array<double> b_;
105 scoped_array<double> D_;
108 #ifndef CERES_NO_SUITESPARSE
109 TEST_F(SparseNormalCholeskySolverTest,
110 SparseNormalCholeskyUsingSuiteSparsePreOrdering) {
111 LinearSolver::Options options;
112 options.sparse_linear_algebra_library_type = SUITE_SPARSE;
113 options.type = SPARSE_NORMAL_CHOLESKY;
114 options.use_postordering = false;
118 TEST_F(SparseNormalCholeskySolverTest,
119 SparseNormalCholeskyUsingSuiteSparsePostOrdering) {
120 LinearSolver::Options options;
121 options.sparse_linear_algebra_library_type = SUITE_SPARSE;
122 options.type = SPARSE_NORMAL_CHOLESKY;
123 options.use_postordering = true;
128 #ifndef CERES_NO_CXSPARSE
129 TEST_F(SparseNormalCholeskySolverTest,
130 SparseNormalCholeskyUsingCXSparsePreOrdering) {
131 LinearSolver::Options options;
132 options.sparse_linear_algebra_library_type = CX_SPARSE;
133 options.type = SPARSE_NORMAL_CHOLESKY;
134 options.use_postordering = false;
138 TEST_F(SparseNormalCholeskySolverTest,
139 SparseNormalCholeskyUsingCXSparsePostOrdering) {
140 LinearSolver::Options options;
141 options.sparse_linear_algebra_library_type = CX_SPARSE;
142 options.type = SPARSE_NORMAL_CHOLESKY;
143 options.use_postordering = true;
148 #ifdef CERES_USE_EIGEN_SPARSE
149 TEST_F(SparseNormalCholeskySolverTest,
150 SparseNormalCholeskyUsingEigenPreOrdering) {
151 LinearSolver::Options options;
152 options.sparse_linear_algebra_library_type = EIGEN_SPARSE;
153 options.type = SPARSE_NORMAL_CHOLESKY;
154 options.use_postordering = false;
158 TEST_F(SparseNormalCholeskySolverTest,
159 SparseNormalCholeskyUsingEigenPostOrdering) {
160 LinearSolver::Options options;
161 options.sparse_linear_algebra_library_type = EIGEN_SPARSE;
162 options.type = SPARSE_NORMAL_CHOLESKY;
163 options.use_postordering = true;
166 #endif // CERES_USE_EIGEN_SPARSE
168 } // namespace internal