Spectra 1.1.0
Header-only C++ Library for Large Scale Eigenvalue Problems
Loading...
Searching...
No Matches
SparseRegularInverse.h
1// Copyright (C) 2017-2025 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_REGULAR_INVERSE_H
8#define SPECTRA_SPARSE_REGULAR_INVERSE_H
9
10#include <Eigen/Core>
11#include <Eigen/SparseCore>
12#include <Eigen/IterativeLinearSolvers>
13#include <stdexcept>
14
15namespace Spectra {
16
36template <typename Scalar_, int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor, typename StorageIndex = int>
38{
39public:
43 using Scalar = Scalar_;
44
45private:
46 using Index = Eigen::Index;
47 using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
48 using MapConstVec = Eigen::Map<const Vector>;
49 using MapVec = Eigen::Map<Vector>;
50 using SparseMatrix = Eigen::SparseMatrix<Scalar, Flags, StorageIndex>;
51 using ConstGenericSparseMatrix = const Eigen::Ref<const SparseMatrix>;
52
53 ConstGenericSparseMatrix m_mat;
54 const Index m_n;
55 Eigen::ConjugateGradient<SparseMatrix> m_cg;
56 mutable CompInfo m_info;
57
58public:
66 template <typename Derived>
67 SparseRegularInverse(const Eigen::SparseMatrixBase<Derived>& mat) :
68 m_mat(mat), m_n(mat.rows())
69 {
70 static_assert(
71 static_cast<int>(Derived::PlainObject::IsRowMajor) == static_cast<int>(SparseMatrix::IsRowMajor),
72 "SparseRegularInverse: 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("SparseRegularInverse: matrix must be square");
76
77 m_cg.compute(mat);
78 m_info = (m_cg.info() == Eigen::Success) ?
81 }
82
86 Index rows() const { return m_n; }
90 Index cols() const { return m_n; }
91
96 CompInfo info() const { return m_info; }
97
104 // y_out = inv(B) * x_in
105 void solve(const Scalar* x_in, Scalar* y_out) const
106 {
107 MapConstVec x(x_in, m_n);
108 MapVec y(y_out, m_n);
109 y.noalias() = m_cg.solve(x);
110
111 m_info = (m_cg.info() == Eigen::Success) ?
114 if (m_info != CompInfo::Successful)
115 throw std::runtime_error("SparseRegularInverse: CG solver does not converge");
116 }
117
124 // y_out = B * x_in
125 void perform_op(const Scalar* x_in, Scalar* y_out) const
126 {
127 MapConstVec x(x_in, m_n);
128 MapVec y(y_out, m_n);
129 y.noalias() = m_mat.template selfadjointView<Uplo>() * x;
130 }
131};
132
133} // namespace Spectra
134
135#endif // SPECTRA_SPARSE_REGULAR_INVERSE_H
void solve(const Scalar *x_in, Scalar *y_out) const
SparseRegularInverse(const Eigen::SparseMatrixBase< Derived > &mat)
void perform_op(const Scalar *x_in, Scalar *y_out) const
@ Successful
Computation was successful.
Definition CompInfo.h:19