Spectra  1.0.1
Header-only C++ Library for Large Scale Eigenvalue Problems
SparseGenRealShiftSolve.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_SPARSE_GEN_REAL_SHIFT_SOLVE_H
8 #define SPECTRA_SPARSE_GEN_REAL_SHIFT_SOLVE_H
9 
10 #include <Eigen/Core>
11 #include <Eigen/SparseCore>
12 #include <Eigen/SparseLU>
13 #include <stdexcept>
14 
15 namespace Spectra {
16 
30 template <typename Scalar_, int Flags = Eigen::ColMajor, typename StorageIndex = int>
32 {
33 public:
37  using Scalar = Scalar_;
38 
39 private:
40  using Index = Eigen::Index;
41  using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
42  using MapConstVec = Eigen::Map<const Vector>;
43  using MapVec = Eigen::Map<Vector>;
44  using SparseMatrix = Eigen::SparseMatrix<Scalar, Flags, StorageIndex>;
45  using ConstGenericSparseMatrix = const Eigen::Ref<const SparseMatrix>;
46 
47  ConstGenericSparseMatrix m_mat;
48  const Index m_n;
49  Eigen::SparseLU<SparseMatrix> m_solver;
50 
51 public:
59  template <typename Derived>
60  SparseGenRealShiftSolve(const Eigen::SparseMatrixBase<Derived>& mat) :
61  m_mat(mat), m_n(mat.rows())
62  {
63  static_assert(
64  static_cast<int>(Derived::PlainObject::IsRowMajor) == static_cast<int>(SparseMatrix::IsRowMajor),
65  "SparseGenRealShiftSolve: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
66 
67  if (mat.rows() != mat.cols())
68  throw std::invalid_argument("SparseGenRealShiftSolve: matrix must be square");
69  }
70 
74  Index rows() const { return m_n; }
78  Index cols() const { return m_n; }
79 
83  void set_shift(const Scalar& sigma)
84  {
85  SparseMatrix I(m_n, m_n);
86  I.setIdentity();
87 
88  m_solver.compute(m_mat - sigma * I);
89  if (m_solver.info() != Eigen::Success)
90  throw std::invalid_argument("SparseGenRealShiftSolve: factorization failed with the given shift");
91  }
92 
99  // y_out = inv(A - sigma * I) * x_in
100  void perform_op(const Scalar* x_in, Scalar* y_out) const
101  {
102  MapConstVec x(x_in, m_n);
103  MapVec y(y_out, m_n);
104  y.noalias() = m_solver.solve(x);
105  }
106 };
107 
108 } // namespace Spectra
109 
110 #endif // SPECTRA_SPARSE_GEN_REAL_SHIFT_SOLVE_H
void perform_op(const Scalar *x_in, Scalar *y_out) const
SparseGenRealShiftSolve(const Eigen::SparseMatrixBase< Derived > &mat)