Spectra 1.1.0
Header-only C++ Library for Large Scale Eigenvalue Problems
Loading...
Searching...
No Matches
SymShiftInvert.h
1// Copyright (C) 2020-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_SYM_SHIFT_INVERT_H
8#define SPECTRA_SYM_SHIFT_INVERT_H
9
10#include <Eigen/Core>
11#include <Eigen/SparseCore>
12#include <Eigen/SparseLU>
13#include <stdexcept>
14#include <type_traits> // std::conditional, std::is_same
15
16#include "../LinAlg/BKLDLT.h"
17#include "../Util/CompInfo.h"
18
19namespace Spectra {
20
22
23// Compute and factorize A-sigma*B without unnecessary copying
24// Default case: A is sparse, B is sparse
25template <bool AIsSparse, bool BIsSparse, int UploA, int UploB>
26class SymShiftInvertHelper
27{
28public:
29 template <typename Scalar, typename Fac, typename ArgA, typename ArgB>
30 static bool factorize(Fac& fac, const ArgA& A, const ArgB& B, const Scalar& sigma)
31 {
32 using SpMat = typename ArgA::PlainObject;
33 SpMat matA = A.template selfadjointView<UploA>();
34 SpMat matB = B.template selfadjointView<UploB>();
35 SpMat mat = matA - sigma * matB;
36 // SparseLU solver
37 fac.isSymmetric(true);
38 fac.compute(mat);
39 // Return true if successful
40 return fac.info() == Eigen::Success;
41 }
42};
43
44// A is dense, B is dense or sparse
45template <bool BIsSparse, int UploA, int UploB>
46class SymShiftInvertHelper<false, BIsSparse, UploA, UploB>
47{
48public:
49 template <typename Scalar, typename Fac, typename ArgA, typename ArgB>
50 static bool factorize(Fac& fac, const ArgA& A, const ArgB& B, const Scalar& sigma)
51 {
52 using Matrix = typename ArgA::PlainObject;
53 // Make a copy of the <UploA> triangular part of A
54 Matrix mat(A.rows(), A.cols());
55 mat.template triangularView<UploA>() = A;
56 // Update <UploA> triangular part of mat
57 if (UploA == UploB)
58 mat -= (B * sigma).template triangularView<UploA>();
59 else
60 mat -= (B * sigma).template triangularView<UploB>().transpose();
61 // BKLDLT solver
62 fac.compute(mat, UploA);
63 // Return true if successful
64 return fac.info() == CompInfo::Successful;
65 }
66};
67
68// A is sparse, B is dense
69template <int UploA, int UploB>
70class SymShiftInvertHelper<true, false, UploA, UploB>
71{
72public:
73 template <typename Scalar, typename Fac, typename ArgA, typename ArgB>
74 static bool factorize(Fac& fac, const ArgA& A, const ArgB& B, const Scalar& sigma)
75 {
76 using Matrix = typename ArgB::PlainObject;
77 // Construct the <UploB> triangular part of -sigma*B
78 Matrix mat(B.rows(), B.cols());
79 mat.template triangularView<UploB>() = -sigma * B;
80 // Update <UploB> triangular part of mat
81 if (UploA == UploB)
82 mat += A.template triangularView<UploB>();
83 else
84 mat += A.template triangularView<UploA>().transpose();
85 // BKLDLT solver
86 fac.compute(mat, UploB);
87 // Return true if successful
88 return fac.info() == CompInfo::Successful;
89 }
90};
91
93
123template <typename Scalar_, typename TypeA = Eigen::Sparse, typename TypeB = Eigen::Sparse,
124 int UploA = Eigen::Lower, int UploB = Eigen::Lower,
125 int FlagsA = Eigen::ColMajor, int FlagsB = Eigen::ColMajor,
126 typename StorageIndexA = int, typename StorageIndexB = int>
128{
129public:
133 using Scalar = Scalar_;
134
135private:
136 using Index = Eigen::Index;
137
138 // Hypothetical type of the A matrix, either dense or sparse
139 using DenseTypeA = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, FlagsA>;
140 using SparseTypeA = Eigen::SparseMatrix<Scalar, FlagsA, StorageIndexA>;
141 // Whether A is sparse
142 using ASparse = std::is_same<TypeA, Eigen::Sparse>;
143 // Actual type of the A matrix
144 using MatrixA = typename std::conditional<ASparse::value, SparseTypeA, DenseTypeA>::type;
145
146 // Hypothetical type of the B matrix, either dense or sparse
147 using DenseTypeB = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, FlagsB>;
148 using SparseTypeB = Eigen::SparseMatrix<Scalar, FlagsB, StorageIndexB>;
149 // Whether B is sparse
150 using BSparse = std::is_same<TypeB, Eigen::Sparse>;
151 // Actual type of the B matrix
152 using MatrixB = typename std::conditional<BSparse::value, SparseTypeB, DenseTypeB>::type;
153
154 using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
155 using MapConstVec = Eigen::Map<const Vector>;
156 using MapVec = Eigen::Map<Vector>;
157
158 // The type of A-sigma*B if one of A and B is dense
159 // DenseType = if (A is dense) MatrixA else MatrixB
160 using DenseType = typename std::conditional<ASparse::value, MatrixB, MatrixA>::type;
161 // The type of A-sigma*B
162 // If both A and B are sparse, the result is MatrixA; otherwise the result is DenseType
163 using ResType = typename std::conditional<ASparse::value && BSparse::value, MatrixA, DenseType>::type;
164
165 // If both A and B are sparse, then the result A-sigma*B is sparse, so we use
166 // sparseLU for factorization; otherwise A-sigma*B is dense, and we use BKLDLT
167 using FacType = typename std::conditional<
168 ASparse::value && BSparse::value,
169 Eigen::SparseLU<ResType>,
170 BKLDLT<Scalar>>::type;
171
172 using ConstGenericMatrixA = const Eigen::Ref<const MatrixA>;
173 using ConstGenericMatrixB = const Eigen::Ref<const MatrixB>;
174
175 ConstGenericMatrixA m_matA;
176 ConstGenericMatrixB m_matB;
177 const Index m_n;
178 FacType m_solver;
179
180public:
190 template <typename DerivedA, typename DerivedB>
191 SymShiftInvert(const Eigen::EigenBase<DerivedA>& A, const Eigen::EigenBase<DerivedB>& B) :
192 m_matA(A.derived()), m_matB(B.derived()), m_n(A.rows())
193 {
194 static_assert(
195 static_cast<int>(DerivedA::PlainObject::IsRowMajor) == static_cast<int>(MatrixA::IsRowMajor),
196 "SymShiftInvert: the \"FlagsA\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
197
198 static_assert(
199 static_cast<int>(DerivedB::PlainObject::IsRowMajor) == static_cast<int>(MatrixB::IsRowMajor),
200 "SymShiftInvert: the \"FlagsB\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
201
202 if (m_n != A.cols() || m_n != B.rows() || m_n != B.cols())
203 throw std::invalid_argument("SymShiftInvert: A and B must be square matrices of the same size");
204 }
205
209 Index rows() const { return m_n; }
213 Index cols() const { return m_n; }
214
218 void set_shift(const Scalar& sigma)
219 {
220 constexpr bool AIsSparse = ASparse::value;
221 constexpr bool BIsSparse = BSparse::value;
222 using Helper = SymShiftInvertHelper<AIsSparse, BIsSparse, UploA, UploB>;
223 const bool success = Helper::factorize(m_solver, m_matA, m_matB, sigma);
224 if (!success)
225 throw std::invalid_argument("SymShiftInvert: factorization failed with the given shift");
226 }
227
234 // y_out = inv(A - sigma * B) * x_in
235 void perform_op(const Scalar* x_in, Scalar* y_out) const
236 {
237 MapConstVec x(x_in, m_n);
238 MapVec y(y_out, m_n);
239 y.noalias() = m_solver.solve(x);
240 }
241};
242
243} // namespace Spectra
244
245#endif // SPECTRA_SYM_SHIFT_INVERT_H
SymShiftInvert(const Eigen::EigenBase< DerivedA > &A, const Eigen::EigenBase< DerivedB > &B)
void set_shift(const Scalar &sigma)
void perform_op(const Scalar *x_in, Scalar *y_out) const