dense_cholesky.h 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308
  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. #ifndef CERES_INTERNAL_DENSE_CHOLESKY_H_
  31. #define CERES_INTERNAL_DENSE_CHOLESKY_H_
  32. // This include must come before any #ifndef check on Ceres compile options.
  33. // clang-format off
  34. #include "ceres/internal/config.h"
  35. // clang-format on
  36. #include <memory>
  37. #include <vector>
  38. #include "Eigen/Dense"
  39. #include "ceres/context_impl.h"
  40. #include "ceres/cuda_buffer.h"
  41. #include "ceres/linear_solver.h"
  42. #include "glog/logging.h"
  43. #ifndef CERES_NO_CUDA
  44. #include "ceres/context_impl.h"
  45. #include "cuda_runtime.h"
  46. #include "cusolverDn.h"
  47. #endif // CERES_NO_CUDA
  48. namespace ceres::internal {
  49. // An interface that abstracts away the internal details of various dense linear
  50. // algebra libraries and offers a simple API for solving dense symmetric
  51. // positive definite linear systems using a Cholesky factorization.
  52. class CERES_NO_EXPORT DenseCholesky {
  53. public:
  54. static std::unique_ptr<DenseCholesky> Create(
  55. const LinearSolver::Options& options);
  56. virtual ~DenseCholesky();
  57. // Computes the Cholesky factorization of the given matrix.
  58. //
  59. // The input matrix lhs is assumed to be a column-major num_cols x num_cols
  60. // matrix, that is symmetric positive definite with its lower triangular part
  61. // containing the left hand side of the linear system being solved.
  62. //
  63. // The input matrix lhs may be modified by the implementation to store the
  64. // factorization, irrespective of whether the factorization succeeds or not.
  65. // As a result it is the user's responsibility to ensure that lhs is valid
  66. // when Solve is called.
  67. virtual LinearSolverTerminationType Factorize(int num_cols,
  68. double* lhs,
  69. std::string* message) = 0;
  70. // Computes the solution to the equation
  71. //
  72. // lhs * solution = rhs
  73. //
  74. // Calling Solve without calling Factorize is undefined behaviour. It is the
  75. // user's responsibility to ensure that the input matrix lhs passed to
  76. // Factorize has not been freed/modified when Solve is called.
  77. virtual LinearSolverTerminationType Solve(const double* rhs,
  78. double* solution,
  79. std::string* message) = 0;
  80. // Convenience method which combines a call to Factorize and Solve. Solve is
  81. // only called if Factorize returns LinearSolverTerminationType::SUCCESS.
  82. //
  83. // The input matrix lhs may be modified by the implementation to store the
  84. // factorization, irrespective of whether the method succeeds or not. It is
  85. // the user's responsibility to ensure that lhs is valid if and when Solve is
  86. // called again after this call.
  87. LinearSolverTerminationType FactorAndSolve(int num_cols,
  88. double* lhs,
  89. const double* rhs,
  90. double* solution,
  91. std::string* message);
  92. };
  93. class CERES_NO_EXPORT EigenDenseCholesky final : public DenseCholesky {
  94. public:
  95. LinearSolverTerminationType Factorize(int num_cols,
  96. double* lhs,
  97. std::string* message) override;
  98. LinearSolverTerminationType Solve(const double* rhs,
  99. double* solution,
  100. std::string* message) override;
  101. private:
  102. using LLTType = Eigen::LLT<Eigen::Ref<Eigen::MatrixXd>, Eigen::Lower>;
  103. std::unique_ptr<LLTType> llt_;
  104. };
  105. class CERES_NO_EXPORT FloatEigenDenseCholesky final : public DenseCholesky {
  106. public:
  107. LinearSolverTerminationType Factorize(int num_cols,
  108. double* lhs,
  109. std::string* message) override;
  110. LinearSolverTerminationType Solve(const double* rhs,
  111. double* solution,
  112. std::string* message) override;
  113. private:
  114. Eigen::MatrixXf lhs_;
  115. Eigen::VectorXf rhs_;
  116. Eigen::VectorXf solution_;
  117. using LLTType = Eigen::LLT<Eigen::MatrixXf, Eigen::Lower>;
  118. std::unique_ptr<LLTType> llt_;
  119. };
  120. #ifndef CERES_NO_LAPACK
  121. class CERES_NO_EXPORT LAPACKDenseCholesky final : public DenseCholesky {
  122. public:
  123. LinearSolverTerminationType Factorize(int num_cols,
  124. double* lhs,
  125. std::string* message) override;
  126. LinearSolverTerminationType Solve(const double* rhs,
  127. double* solution,
  128. std::string* message) override;
  129. private:
  130. double* lhs_ = nullptr;
  131. int num_cols_ = -1;
  132. LinearSolverTerminationType termination_type_ =
  133. LinearSolverTerminationType::FATAL_ERROR;
  134. };
  135. class CERES_NO_EXPORT FloatLAPACKDenseCholesky final : public DenseCholesky {
  136. public:
  137. LinearSolverTerminationType Factorize(int num_cols,
  138. double* lhs,
  139. std::string* message) override;
  140. LinearSolverTerminationType Solve(const double* rhs,
  141. double* solution,
  142. std::string* message) override;
  143. private:
  144. Eigen::MatrixXf lhs_;
  145. Eigen::VectorXf rhs_and_solution_;
  146. int num_cols_ = -1;
  147. LinearSolverTerminationType termination_type_ =
  148. LinearSolverTerminationType::FATAL_ERROR;
  149. };
  150. #endif // CERES_NO_LAPACK
  151. class DenseIterativeRefiner;
  152. // Computes an initial solution using the given instance of
  153. // DenseCholesky, and then refines it using the DenseIterativeRefiner.
  154. class CERES_NO_EXPORT RefinedDenseCholesky final : public DenseCholesky {
  155. public:
  156. RefinedDenseCholesky(
  157. std::unique_ptr<DenseCholesky> dense_cholesky,
  158. std::unique_ptr<DenseIterativeRefiner> iterative_refiner);
  159. ~RefinedDenseCholesky() override;
  160. LinearSolverTerminationType Factorize(int num_cols,
  161. double* lhs,
  162. std::string* message) override;
  163. LinearSolverTerminationType Solve(const double* rhs,
  164. double* solution,
  165. std::string* message) override;
  166. private:
  167. std::unique_ptr<DenseCholesky> dense_cholesky_;
  168. std::unique_ptr<DenseIterativeRefiner> iterative_refiner_;
  169. double* lhs_ = nullptr;
  170. int num_cols_;
  171. };
  172. #ifndef CERES_NO_CUDA
  173. // CUDA implementation of DenseCholesky using the cuSolverDN library using the
  174. // 32-bit legacy interface for maximum compatibility.
  175. class CERES_NO_EXPORT CUDADenseCholesky final : public DenseCholesky {
  176. public:
  177. static std::unique_ptr<CUDADenseCholesky> Create(
  178. const LinearSolver::Options& options);
  179. CUDADenseCholesky(const CUDADenseCholesky&) = delete;
  180. CUDADenseCholesky& operator=(const CUDADenseCholesky&) = delete;
  181. LinearSolverTerminationType Factorize(int num_cols,
  182. double* lhs,
  183. std::string* message) override;
  184. LinearSolverTerminationType Solve(const double* rhs,
  185. double* solution,
  186. std::string* message) override;
  187. private:
  188. explicit CUDADenseCholesky(ContextImpl* context);
  189. ContextImpl* context_ = nullptr;
  190. // Number of columns in the A matrix, to be cached between calls to *Factorize
  191. // and *Solve.
  192. size_t num_cols_ = 0;
  193. // GPU memory allocated for the A matrix (lhs matrix).
  194. CudaBuffer<double> lhs_;
  195. // GPU memory allocated for the B matrix (rhs vector).
  196. CudaBuffer<double> rhs_;
  197. // Scratch space for cuSOLVER on the GPU.
  198. CudaBuffer<double> device_workspace_;
  199. // Required for error handling with cuSOLVER.
  200. CudaBuffer<int> error_;
  201. // Cache the result of Factorize to ensure that when Solve is called, the
  202. // factorization of lhs is valid.
  203. LinearSolverTerminationType factorize_result_ =
  204. LinearSolverTerminationType::FATAL_ERROR;
  205. };
  206. // A mixed-precision iterative refinement dense Cholesky solver using FP32 CUDA
  207. // Dense Cholesky for inner iterations, and FP64 outer refinements.
  208. // This class implements a modified version of the "Classical iterative
  209. // refinement" (Algorithm 4.1) from the following paper:
  210. // Haidar, Azzam, Harun Bayraktar, Stanimire Tomov, Jack Dongarra, and Nicholas
  211. // J. Higham. "Mixed-precision iterative refinement using tensor cores on GPUs
  212. // to accelerate solution of linear systems." Proceedings of the Royal Society A
  213. // 476, no. 2243 (2020): 20200110.
  214. //
  215. // The three key modifications from Algorithm 4.1 in the paper are:
  216. // 1. We use Cholesky factorization instead of LU factorization since our A is
  217. // symmetric positive definite.
  218. // 2. During the solution update, the up-cast and accumulation is performed in
  219. // one step with a custom kernel.
  220. class CERES_NO_EXPORT CUDADenseCholeskyMixedPrecision final
  221. : public DenseCholesky {
  222. public:
  223. static std::unique_ptr<CUDADenseCholeskyMixedPrecision> Create(
  224. const LinearSolver::Options& options);
  225. CUDADenseCholeskyMixedPrecision(const CUDADenseCholeskyMixedPrecision&) =
  226. delete;
  227. CUDADenseCholeskyMixedPrecision& operator=(
  228. const CUDADenseCholeskyMixedPrecision&) = delete;
  229. LinearSolverTerminationType Factorize(int num_cols,
  230. double* lhs,
  231. std::string* message) override;
  232. LinearSolverTerminationType Solve(const double* rhs,
  233. double* solution,
  234. std::string* message) override;
  235. private:
  236. CUDADenseCholeskyMixedPrecision(ContextImpl* context,
  237. int max_num_refinement_iterations);
  238. // Helper function to wrap Cuda boilerplate needed to call Spotrf.
  239. LinearSolverTerminationType CudaCholeskyFactorize(std::string* message);
  240. // Helper function to wrap Cuda boilerplate needed to call Spotrs.
  241. LinearSolverTerminationType CudaCholeskySolve(std::string* message);
  242. // Picks up the cuSolverDN and cuStream handles from the context in the
  243. // options, and the number of refinement iterations from the options. If
  244. // the context is unable to initialize CUDA, returns false with a
  245. // human-readable message indicating the reason.
  246. bool Init(const LinearSolver::Options& options, std::string* message);
  247. ContextImpl* context_ = nullptr;
  248. // Number of columns in the A matrix, to be cached between calls to *Factorize
  249. // and *Solve.
  250. size_t num_cols_ = 0;
  251. CudaBuffer<double> lhs_fp64_;
  252. CudaBuffer<double> rhs_fp64_;
  253. CudaBuffer<float> lhs_fp32_;
  254. // Scratch space for cuSOLVER on the GPU.
  255. CudaBuffer<float> device_workspace_;
  256. // Required for error handling with cuSOLVER.
  257. CudaBuffer<int> error_;
  258. // Solution to lhs * x = rhs.
  259. CudaBuffer<double> x_fp64_;
  260. // Incremental correction to x.
  261. CudaBuffer<float> correction_fp32_;
  262. // Residual to iterative refinement.
  263. CudaBuffer<float> residual_fp32_;
  264. CudaBuffer<double> residual_fp64_;
  265. // Number of inner refinement iterations to perform.
  266. int max_num_refinement_iterations_ = 0;
  267. // Cache the result of Factorize to ensure that when Solve is called, the
  268. // factorization of lhs is valid.
  269. LinearSolverTerminationType factorize_result_ =
  270. LinearSolverTerminationType::FATAL_ERROR;
  271. };
  272. #endif // CERES_NO_CUDA
  273. } // namespace ceres::internal
  274. #endif // CERES_INTERNAL_DENSE_CHOLESKY_H_