1 // Ceres Solver - A fast non-linear least squares minimizer 2 // Copyright 2010, 2011, 2012 Google Inc. All rights reserved. 3 // http://code.google.com/p/ceres-solver/ 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 (at) google.com (Sameer Agarwal) 30 31 #include "ceres/casts.h" 32 #include "ceres/compressed_row_sparse_matrix.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" 40 41 42 namespace ceres { 43 namespace internal { 44 45 class UnsymmetricLinearSolverTest : public ::testing::Test { 46 protected : 47 virtual void SetUp() { 48 scoped_ptr<LinearLeastSquaresProblem> problem( 49 CreateLinearLeastSquaresProblemFromId(0)); 50 51 CHECK_NOTNULL(problem.get()); 52 A_.reset(down_cast<TripletSparseMatrix*>(problem->A.release())); 53 b_.reset(problem->b.release()); 54 D_.reset(problem->D.release()); 55 sol_unregularized_.reset(problem->x.release()); 56 sol_regularized_.reset(problem->x_D.release()); 57 } 58 59 void TestSolver(const LinearSolver::Options& options) { 60 scoped_ptr<LinearSolver> solver(LinearSolver::Create(options)); 61 62 LinearSolver::PerSolveOptions per_solve_options; 63 LinearSolver::Summary unregularized_solve_summary; 64 LinearSolver::Summary regularized_solve_summary; 65 Vector x_unregularized(A_->num_cols()); 66 Vector x_regularized(A_->num_cols()); 67 68 scoped_ptr<SparseMatrix> transformed_A; 69 70 if (options.type == DENSE_QR || 71 options.type == DENSE_NORMAL_CHOLESKY) { 72 transformed_A.reset(new DenseSparseMatrix(*A_)); 73 } else if (options.type == SPARSE_NORMAL_CHOLESKY) { 74 CompressedRowSparseMatrix* crsm = new CompressedRowSparseMatrix(*A_); 75 // Add row/column blocks structure. 76 for (int i = 0; i < A_->num_rows(); ++i) { 77 crsm->mutable_row_blocks()->push_back(1); 78 } 79 80 for (int i = 0; i < A_->num_cols(); ++i) { 81 crsm->mutable_col_blocks()->push_back(1); 82 } 83 transformed_A.reset(crsm); 84 } else { 85 LOG(FATAL) << "Unknown linear solver : " << options.type; 86 } 87 // Unregularized 88 unregularized_solve_summary = 89 solver->Solve(transformed_A.get(), 90 b_.get(), 91 per_solve_options, 92 x_unregularized.data()); 93 94 // Regularized solution 95 per_solve_options.D = D_.get(); 96 regularized_solve_summary = 97 solver->Solve(transformed_A.get(), 98 b_.get(), 99 per_solve_options, 100 x_regularized.data()); 101 102 EXPECT_EQ(unregularized_solve_summary.termination_type, TOLERANCE); 103 104 for (int i = 0; i < A_->num_cols(); ++i) { 105 EXPECT_NEAR(sol_unregularized_[i], x_unregularized[i], 1e-8); 106 } 107 108 EXPECT_EQ(regularized_solve_summary.termination_type, TOLERANCE); 109 for (int i = 0; i < A_->num_cols(); ++i) { 110 EXPECT_NEAR(sol_regularized_[i], x_regularized[i], 1e-8); 111 } 112 } 113 114 scoped_ptr<TripletSparseMatrix> A_; 115 scoped_array<double> b_; 116 scoped_array<double> D_; 117 scoped_array<double> sol_unregularized_; 118 scoped_array<double> sol_regularized_; 119 }; 120 121 TEST_F(UnsymmetricLinearSolverTest, EigenDenseQR) { 122 LinearSolver::Options options; 123 options.type = DENSE_QR; 124 options.dense_linear_algebra_library_type = EIGEN; 125 TestSolver(options); 126 } 127 128 TEST_F(UnsymmetricLinearSolverTest, EigenDenseNormalCholesky) { 129 LinearSolver::Options options; 130 options.dense_linear_algebra_library_type = EIGEN; 131 options.type = DENSE_NORMAL_CHOLESKY; 132 TestSolver(options); 133 } 134 135 #ifndef CERES_NO_LAPACK 136 TEST_F(UnsymmetricLinearSolverTest, LAPACKDenseQR) { 137 LinearSolver::Options options; 138 options.type = DENSE_QR; 139 options.dense_linear_algebra_library_type = LAPACK; 140 TestSolver(options); 141 } 142 143 TEST_F(UnsymmetricLinearSolverTest, LAPACKDenseNormalCholesky) { 144 LinearSolver::Options options; 145 options.dense_linear_algebra_library_type = LAPACK; 146 options.type = DENSE_NORMAL_CHOLESKY; 147 TestSolver(options); 148 } 149 #endif 150 151 #ifndef CERES_NO_SUITESPARSE 152 TEST_F(UnsymmetricLinearSolverTest, 153 SparseNormalCholeskyUsingSuiteSparsePreOrdering) { 154 LinearSolver::Options options; 155 options.sparse_linear_algebra_library_type = SUITE_SPARSE; 156 options.type = SPARSE_NORMAL_CHOLESKY; 157 options.use_postordering = false; 158 TestSolver(options); 159 } 160 161 TEST_F(UnsymmetricLinearSolverTest, 162 SparseNormalCholeskyUsingSuiteSparsePostOrdering) { 163 LinearSolver::Options options; 164 options.sparse_linear_algebra_library_type = SUITE_SPARSE; 165 options.type = SPARSE_NORMAL_CHOLESKY; 166 options.use_postordering = true; 167 TestSolver(options); 168 } 169 #endif 170 171 #ifndef CERES_NO_CXSPARSE 172 TEST_F(UnsymmetricLinearSolverTest, 173 SparseNormalCholeskyUsingCXSparsePreOrdering) { 174 LinearSolver::Options options; 175 options.sparse_linear_algebra_library_type = CX_SPARSE; 176 options.type = SPARSE_NORMAL_CHOLESKY; 177 options.use_postordering = false; 178 TestSolver(options); 179 } 180 181 TEST_F(UnsymmetricLinearSolverTest, 182 SparseNormalCholeskyUsingCXSparsePostOrdering) { 183 LinearSolver::Options options; 184 options.sparse_linear_algebra_library_type = CX_SPARSE; 185 options.type = SPARSE_NORMAL_CHOLESKY; 186 options.use_postordering = true; 187 TestSolver(options); 188 } 189 #endif 190 191 } // namespace internal 192 } // namespace ceres 193