Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 1 | // 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/eigensparse.h" |
| 32 | |
| 33 | #ifdef CERES_USE_EIGEN_SPARSE |
| 34 | |
| 35 | #include <sstream> |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 36 | |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 37 | #include "Eigen/SparseCholesky" |
| 38 | #include "Eigen/SparseCore" |
| 39 | #include "ceres/compressed_row_sparse_matrix.h" |
| 40 | #include "ceres/linear_solver.h" |
| 41 | |
| 42 | namespace ceres { |
| 43 | namespace internal { |
| 44 | |
| 45 | // TODO(sameeragarwal): Use enable_if to clean up the implementations |
| 46 | // for when Scalar == double. |
| 47 | template <typename Solver> |
| 48 | class EigenSparseCholeskyTemplate : public SparseCholesky { |
| 49 | public: |
| 50 | EigenSparseCholeskyTemplate() : analyzed_(false) {} |
| 51 | virtual ~EigenSparseCholeskyTemplate() {} |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 52 | CompressedRowSparseMatrix::StorageType StorageType() const final { |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 53 | return CompressedRowSparseMatrix::LOWER_TRIANGULAR; |
| 54 | } |
| 55 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 56 | LinearSolverTerminationType Factorize( |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 57 | const Eigen::SparseMatrix<typename Solver::Scalar>& lhs, |
| 58 | std::string* message) { |
| 59 | if (!analyzed_) { |
| 60 | solver_.analyzePattern(lhs); |
| 61 | |
| 62 | if (VLOG_IS_ON(2)) { |
| 63 | std::stringstream ss; |
| 64 | solver_.dumpMemory(ss); |
| 65 | VLOG(2) << "Symbolic Analysis\n" << ss.str(); |
| 66 | } |
| 67 | |
| 68 | if (solver_.info() != Eigen::Success) { |
| 69 | *message = "Eigen failure. Unable to find symbolic factorization."; |
| 70 | return LINEAR_SOLVER_FATAL_ERROR; |
| 71 | } |
| 72 | |
| 73 | analyzed_ = true; |
| 74 | } |
| 75 | |
| 76 | solver_.factorize(lhs); |
| 77 | if (solver_.info() != Eigen::Success) { |
| 78 | *message = "Eigen failure. Unable to find numeric factorization."; |
| 79 | return LINEAR_SOLVER_FAILURE; |
| 80 | } |
| 81 | return LINEAR_SOLVER_SUCCESS; |
| 82 | } |
| 83 | |
| 84 | LinearSolverTerminationType Solve(const double* rhs_ptr, |
| 85 | double* solution_ptr, |
| 86 | std::string* message) { |
| 87 | CHECK(analyzed_) << "Solve called without a call to Factorize first."; |
| 88 | |
| 89 | scalar_rhs_ = ConstVectorRef(rhs_ptr, solver_.cols()) |
| 90 | .template cast<typename Solver::Scalar>(); |
| 91 | |
| 92 | // The two casts are needed if the Scalar in this class is not |
| 93 | // double. For code simplicity we are going to assume that Eigen |
| 94 | // is smart enough to figure out that casting a double Vector to a |
| 95 | // double Vector is a straight copy. If this turns into a |
| 96 | // performance bottleneck (unlikely), we can revisit this. |
| 97 | scalar_solution_ = solver_.solve(scalar_rhs_); |
| 98 | VectorRef(solution_ptr, solver_.cols()) = |
| 99 | scalar_solution_.template cast<double>(); |
| 100 | |
| 101 | if (solver_.info() != Eigen::Success) { |
| 102 | *message = "Eigen failure. Unable to do triangular solve."; |
| 103 | return LINEAR_SOLVER_FAILURE; |
| 104 | } |
| 105 | return LINEAR_SOLVER_SUCCESS; |
| 106 | } |
| 107 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 108 | LinearSolverTerminationType Factorize(CompressedRowSparseMatrix* lhs, |
| 109 | std::string* message) final { |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 110 | CHECK_EQ(lhs->storage_type(), StorageType()); |
| 111 | |
| 112 | typename Solver::Scalar* values_ptr = NULL; |
| 113 | if (std::is_same<typename Solver::Scalar, double>::value) { |
| 114 | values_ptr = |
| 115 | reinterpret_cast<typename Solver::Scalar*>(lhs->mutable_values()); |
| 116 | } else { |
| 117 | // In the case where the scalar used in this class is not |
| 118 | // double. In that case, make a copy of the values array in the |
| 119 | // CompressedRowSparseMatrix and cast it to Scalar along the way. |
| 120 | values_ = ConstVectorRef(lhs->values(), lhs->num_nonzeros()) |
| 121 | .cast<typename Solver::Scalar>(); |
| 122 | values_ptr = values_.data(); |
| 123 | } |
| 124 | |
| 125 | Eigen::MappedSparseMatrix<typename Solver::Scalar, Eigen::ColMajor> |
| 126 | eigen_lhs(lhs->num_rows(), |
| 127 | lhs->num_rows(), |
| 128 | lhs->num_nonzeros(), |
| 129 | lhs->mutable_rows(), |
| 130 | lhs->mutable_cols(), |
| 131 | values_ptr); |
| 132 | return Factorize(eigen_lhs, message); |
| 133 | } |
| 134 | |
| 135 | private: |
| 136 | Eigen::Matrix<typename Solver::Scalar, Eigen::Dynamic, 1> values_, |
| 137 | scalar_rhs_, scalar_solution_; |
| 138 | bool analyzed_; |
| 139 | Solver solver_; |
| 140 | }; |
| 141 | |
| 142 | std::unique_ptr<SparseCholesky> EigenSparseCholesky::Create( |
| 143 | const OrderingType ordering_type) { |
| 144 | std::unique_ptr<SparseCholesky> sparse_cholesky; |
| 145 | |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 146 | typedef Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>, |
| 147 | Eigen::Upper, |
| 148 | Eigen::AMDOrdering<int>> |
| 149 | WithAMDOrdering; |
| 150 | typedef Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>, |
| 151 | Eigen::Upper, |
| 152 | Eigen::NaturalOrdering<int>> |
| 153 | WithNaturalOrdering; |
| 154 | if (ordering_type == AMD) { |
| 155 | sparse_cholesky.reset(new EigenSparseCholeskyTemplate<WithAMDOrdering>()); |
| 156 | } else { |
| 157 | sparse_cholesky.reset( |
| 158 | new EigenSparseCholeskyTemplate<WithNaturalOrdering>()); |
| 159 | } |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 160 | return sparse_cholesky; |
| 161 | } |
| 162 | |
| 163 | EigenSparseCholesky::~EigenSparseCholesky() {} |
| 164 | |
| 165 | std::unique_ptr<SparseCholesky> FloatEigenSparseCholesky::Create( |
| 166 | const OrderingType ordering_type) { |
| 167 | std::unique_ptr<SparseCholesky> sparse_cholesky; |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 168 | typedef Eigen::SimplicialLDLT<Eigen::SparseMatrix<float>, |
| 169 | Eigen::Upper, |
| 170 | Eigen::AMDOrdering<int>> |
| 171 | WithAMDOrdering; |
| 172 | typedef Eigen::SimplicialLDLT<Eigen::SparseMatrix<float>, |
| 173 | Eigen::Upper, |
| 174 | Eigen::NaturalOrdering<int>> |
| 175 | WithNaturalOrdering; |
| 176 | if (ordering_type == AMD) { |
| 177 | sparse_cholesky.reset(new EigenSparseCholeskyTemplate<WithAMDOrdering>()); |
| 178 | } else { |
| 179 | sparse_cholesky.reset( |
| 180 | new EigenSparseCholeskyTemplate<WithNaturalOrdering>()); |
| 181 | } |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 182 | return sparse_cholesky; |
| 183 | } |
| 184 | |
| 185 | FloatEigenSparseCholesky::~FloatEigenSparseCholesky() {} |
| 186 | |
| 187 | } // namespace internal |
| 188 | } // namespace ceres |
| 189 | |
| 190 | #endif // CERES_USE_EIGEN_SPARSE |