Spectra  1.0.1
Header-only C++ Library for Large Scale Eigenvalue Problems
DenseGenComplexShiftSolve.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_GEN_COMPLEX_SHIFT_SOLVE_H
8 #define SPECTRA_DENSE_GEN_COMPLEX_SHIFT_SOLVE_H
9 
10 #include <Eigen/Core>
11 #include <Eigen/LU>
12 #include <stdexcept>
13 
14 namespace Spectra {
15 
29 template <typename Scalar_, int Flags = Eigen::ColMajor>
31 {
32 public:
36  using Scalar = Scalar_;
37 
38 private:
39  using Index = Eigen::Index;
40  using Matrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Flags>;
41  using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
42  using MapConstVec = Eigen::Map<const Vector>;
43  using MapVec = Eigen::Map<Vector>;
44  using ConstGenericMatrix = const Eigen::Ref<const Matrix>;
45 
46  using Complex = std::complex<Scalar>;
47  using ComplexMatrix = Eigen::Matrix<Complex, Eigen::Dynamic, Eigen::Dynamic, Flags>;
48  using ComplexVector = Eigen::Matrix<Complex, Eigen::Dynamic, 1>;
49 
50  using ComplexSolver = Eigen::PartialPivLU<ComplexMatrix>;
51 
52  ConstGenericMatrix m_mat;
53  const Index m_n;
54  ComplexSolver m_solver;
55  mutable ComplexVector m_x_cache;
56 
57 public:
66  template <typename Derived>
67  DenseGenComplexShiftSolve(const Eigen::MatrixBase<Derived>& mat) :
68  m_mat(mat), m_n(mat.rows())
69  {
70  static_assert(
71  static_cast<int>(Derived::PlainObject::IsRowMajor) == static_cast<int>(Matrix::IsRowMajor),
72  "DenseGenComplexShiftSolve: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
73 
74  if (mat.rows() != mat.cols())
75  throw std::invalid_argument("DenseGenComplexShiftSolve: matrix must be square");
76  }
77 
81  Index rows() const { return m_n; }
85  Index cols() const { return m_n; }
86 
93  void set_shift(const Scalar& sigmar, const Scalar& sigmai)
94  {
95  m_solver.compute(m_mat.template cast<Complex>() - Complex(sigmar, sigmai) * ComplexMatrix::Identity(m_n, m_n));
96  m_x_cache.resize(m_n);
97  m_x_cache.setZero();
98  }
99 
107  // y_out = Re( inv(A - sigma * I) * x_in )
108  void perform_op(const Scalar* x_in, Scalar* y_out) const
109  {
110  m_x_cache.real() = MapConstVec(x_in, m_n);
111  MapVec y(y_out, m_n);
112  y.noalias() = m_solver.solve(m_x_cache).real();
113  }
114 };
115 
116 } // namespace Spectra
117 
118 #endif // SPECTRA_DENSE_GEN_COMPLEX_SHIFT_SOLVE_H
DenseGenComplexShiftSolve(const Eigen::MatrixBase< Derived > &mat)
void perform_op(const Scalar *x_in, Scalar *y_out) const
void set_shift(const Scalar &sigmar, const Scalar &sigmai)