eigensparse.cc 8.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2023 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. #include "ceres/eigensparse.h"
  31. #include <memory>
  32. #ifdef CERES_USE_EIGEN_SPARSE
  33. #include <sstream>
  34. #ifndef CERES_NO_EIGEN_METIS
  35. #include <iostream> // This is needed because MetisSupport depends on iostream.
  36. #include "Eigen/MetisSupport"
  37. #endif
  38. #include "Eigen/SparseCholesky"
  39. #include "Eigen/SparseCore"
  40. #include "ceres/compressed_row_sparse_matrix.h"
  41. #include "ceres/linear_solver.h"
  42. namespace ceres::internal {
  43. template <typename Solver>
  44. class EigenSparseCholeskyTemplate final : public SparseCholesky {
  45. public:
  46. EigenSparseCholeskyTemplate() = default;
  47. CompressedRowSparseMatrix::StorageType StorageType() const final {
  48. return CompressedRowSparseMatrix::StorageType::LOWER_TRIANGULAR;
  49. }
  50. LinearSolverTerminationType Factorize(
  51. const Eigen::SparseMatrix<typename Solver::Scalar>& lhs,
  52. std::string* message) {
  53. if (!analyzed_) {
  54. solver_.analyzePattern(lhs);
  55. if (VLOG_IS_ON(2)) {
  56. std::stringstream ss;
  57. solver_.dumpMemory(ss);
  58. VLOG(2) << "Symbolic Analysis\n" << ss.str();
  59. }
  60. if (solver_.info() != Eigen::Success) {
  61. *message = "Eigen failure. Unable to find symbolic factorization.";
  62. return LinearSolverTerminationType::FATAL_ERROR;
  63. }
  64. analyzed_ = true;
  65. }
  66. solver_.factorize(lhs);
  67. if (solver_.info() != Eigen::Success) {
  68. *message = "Eigen failure. Unable to find numeric factorization.";
  69. return LinearSolverTerminationType::FAILURE;
  70. }
  71. return LinearSolverTerminationType::SUCCESS;
  72. }
  73. LinearSolverTerminationType Solve(const double* rhs_ptr,
  74. double* solution_ptr,
  75. std::string* message) override {
  76. CHECK(analyzed_) << "Solve called without a call to Factorize first.";
  77. // Avoid copying when the scalar type is double
  78. if constexpr (std::is_same_v<typename Solver::Scalar, double>) {
  79. ConstVectorRef scalar_rhs(rhs_ptr, solver_.cols());
  80. VectorRef(solution_ptr, solver_.cols()) = solver_.solve(scalar_rhs);
  81. } else {
  82. auto scalar_rhs = ConstVectorRef(rhs_ptr, solver_.cols())
  83. .template cast<typename Solver::Scalar>();
  84. auto scalar_solution = solver_.solve(scalar_rhs);
  85. VectorRef(solution_ptr, solver_.cols()) =
  86. scalar_solution.template cast<double>();
  87. }
  88. if (solver_.info() != Eigen::Success) {
  89. *message = "Eigen failure. Unable to do triangular solve.";
  90. return LinearSolverTerminationType::FAILURE;
  91. }
  92. return LinearSolverTerminationType::SUCCESS;
  93. }
  94. LinearSolverTerminationType Factorize(CompressedRowSparseMatrix* lhs,
  95. std::string* message) final {
  96. CHECK_EQ(lhs->storage_type(), StorageType());
  97. typename Solver::Scalar* values_ptr = nullptr;
  98. if constexpr (std::is_same_v<typename Solver::Scalar, double>) {
  99. values_ptr = lhs->mutable_values();
  100. } else {
  101. // In the case where the scalar used in this class is not
  102. // double. In that case, make a copy of the values array in the
  103. // CompressedRowSparseMatrix and cast it to Scalar along the way.
  104. values_ = ConstVectorRef(lhs->values(), lhs->num_nonzeros())
  105. .cast<typename Solver::Scalar>();
  106. values_ptr = values_.data();
  107. }
  108. Eigen::Map<
  109. const Eigen::SparseMatrix<typename Solver::Scalar, Eigen::ColMajor>>
  110. eigen_lhs(lhs->num_rows(),
  111. lhs->num_rows(),
  112. lhs->num_nonzeros(),
  113. lhs->rows(),
  114. lhs->cols(),
  115. values_ptr);
  116. return Factorize(eigen_lhs, message);
  117. }
  118. private:
  119. Eigen::Matrix<typename Solver::Scalar, Eigen::Dynamic, 1> values_;
  120. bool analyzed_{false};
  121. Solver solver_;
  122. };
  123. std::unique_ptr<SparseCholesky> EigenSparseCholesky::Create(
  124. const OrderingType ordering_type) {
  125. using WithAMDOrdering = Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>,
  126. Eigen::Upper,
  127. Eigen::AMDOrdering<int>>;
  128. using WithNaturalOrdering =
  129. Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>,
  130. Eigen::Upper,
  131. Eigen::NaturalOrdering<int>>;
  132. if (ordering_type == OrderingType::AMD) {
  133. return std::make_unique<EigenSparseCholeskyTemplate<WithAMDOrdering>>();
  134. } else if (ordering_type == OrderingType::NESDIS) {
  135. #ifndef CERES_NO_EIGEN_METIS
  136. using WithMetisOrdering = Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>,
  137. Eigen::Upper,
  138. Eigen::MetisOrdering<int>>;
  139. return std::make_unique<EigenSparseCholeskyTemplate<WithMetisOrdering>>();
  140. #else
  141. LOG(FATAL)
  142. << "Congratulations you have found a bug in Ceres Solver. Please "
  143. "report it to the Ceres Solver developers.";
  144. return nullptr;
  145. #endif // CERES_NO_EIGEN_METIS
  146. }
  147. return std::make_unique<EigenSparseCholeskyTemplate<WithNaturalOrdering>>();
  148. }
  149. EigenSparseCholesky::~EigenSparseCholesky() = default;
  150. std::unique_ptr<SparseCholesky> FloatEigenSparseCholesky::Create(
  151. const OrderingType ordering_type) {
  152. using WithAMDOrdering = Eigen::SimplicialLDLT<Eigen::SparseMatrix<float>,
  153. Eigen::Upper,
  154. Eigen::AMDOrdering<int>>;
  155. using WithNaturalOrdering =
  156. Eigen::SimplicialLDLT<Eigen::SparseMatrix<float>,
  157. Eigen::Upper,
  158. Eigen::NaturalOrdering<int>>;
  159. if (ordering_type == OrderingType::AMD) {
  160. return std::make_unique<EigenSparseCholeskyTemplate<WithAMDOrdering>>();
  161. } else if (ordering_type == OrderingType::NESDIS) {
  162. #ifndef CERES_NO_EIGEN_METIS
  163. using WithMetisOrdering = Eigen::SimplicialLDLT<Eigen::SparseMatrix<float>,
  164. Eigen::Upper,
  165. Eigen::MetisOrdering<int>>;
  166. return std::make_unique<EigenSparseCholeskyTemplate<WithMetisOrdering>>();
  167. #else
  168. LOG(FATAL)
  169. << "Congratulations you have found a bug in Ceres Solver. Please "
  170. "report it to the Ceres Solver developers.";
  171. return nullptr;
  172. #endif // CERES_NO_EIGEN_METIS
  173. }
  174. return std::make_unique<EigenSparseCholeskyTemplate<WithNaturalOrdering>>();
  175. }
  176. FloatEigenSparseCholesky::~FloatEigenSparseCholesky() = default;
  177. } // namespace ceres::internal
  178. #endif // CERES_USE_EIGEN_SPARSE