Spectra  1.0.1
Header-only C++ Library for Large Scale Eigenvalue Problems
DenseCholesky.h
1 // Copyright (C) 2016-2022 Yixuan Qiu <yixuan.qiu@cos.name>
2 //
3 // This Source Code Form is subject to the terms of the Mozilla
4 // Public License v. 2.0. If a copy of the MPL was not distributed
5 // with this file, You can obtain one at https://mozilla.org/MPL/2.0/.
6 
7 #ifndef SPECTRA_DENSE_CHOLESKY_H
8 #define SPECTRA_DENSE_CHOLESKY_H
9 
10 #include <Eigen/Core>
11 #include <Eigen/Cholesky>
12 #include <stdexcept>
13 
14 #include "../Util/CompInfo.h"
15 
16 namespace Spectra {
17 
33 template <typename Scalar_, int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor>
35 {
36 public:
40  using Scalar = Scalar_;
41 
42 private:
43  using Index = Eigen::Index;
44  using Matrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Flags>;
45  using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
46  using MapConstVec = Eigen::Map<const Vector>;
47  using MapVec = Eigen::Map<Vector>;
48 
49  const Index m_n;
50  Eigen::LLT<Matrix, Uplo> m_decomp;
51  CompInfo m_info; // status of the decomposition
52 
53 public:
62  template <typename Derived>
63  DenseCholesky(const Eigen::MatrixBase<Derived>& mat) :
64  m_n(mat.rows()), m_info(CompInfo::NotComputed)
65  {
66  static_assert(
67  static_cast<int>(Derived::PlainObject::IsRowMajor) == static_cast<int>(Matrix::IsRowMajor),
68  "DenseCholesky: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
69 
70  if (m_n != mat.cols())
71  throw std::invalid_argument("DenseCholesky: matrix must be square");
72 
73  m_decomp.compute(mat);
74  m_info = (m_decomp.info() == Eigen::Success) ?
77  }
78 
82  Index rows() const { return m_n; }
86  Index cols() const { return m_n; }
87 
92  CompInfo info() const { return m_info; }
93 
100  // y_out = inv(L) * x_in
101  void lower_triangular_solve(const Scalar* x_in, Scalar* y_out) const
102  {
103  MapConstVec x(x_in, m_n);
104  MapVec y(y_out, m_n);
105  y.noalias() = m_decomp.matrixL().solve(x);
106  }
107 
114  // y_out = inv(L') * x_in
115  void upper_triangular_solve(const Scalar* x_in, Scalar* y_out) const
116  {
117  MapConstVec x(x_in, m_n);
118  MapVec y(y_out, m_n);
119  y.noalias() = m_decomp.matrixU().solve(x);
120  }
121 };
122 
123 } // namespace Spectra
124 
125 #endif // SPECTRA_DENSE_CHOLESKY_H
CompInfo info() const
Definition: DenseCholesky.h:92
void upper_triangular_solve(const Scalar *x_in, Scalar *y_out) const
DenseCholesky(const Eigen::MatrixBase< Derived > &mat)
Definition: DenseCholesky.h:63
void lower_triangular_solve(const Scalar *x_in, Scalar *y_out) const
CompInfo
Definition: CompInfo.h:18
@ Successful
Computation was successful.