Go to the documentation of this file.
7 #ifndef SPECTRA_SYM_SHIFT_INVERT_H
8 #define SPECTRA_SYM_SHIFT_INVERT_H
11 #include <Eigen/SparseCore>
12 #include <Eigen/SparseLU>
14 #include <type_traits>
16 #include "../LinAlg/BKLDLT.h"
17 #include "../Util/CompInfo.h"
25 template <
bool AIsSparse,
bool BIsSparse,
int UploA,
int UploB>
26 class SymShiftInvertHelper
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)
32 using SpMat =
typename ArgA::PlainObject;
33 SpMat matA =
A.template selfadjointView<UploA>();
34 SpMat matB =
B.template selfadjointView<UploB>();
37 fac.isSymmetric(
true);
45 template <
bool BIsSparse,
int UploA,
int UploB>
46 class SymShiftInvertHelper<false, BIsSparse, UploA, UploB>
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)
52 using Matrix =
typename ArgA::PlainObject;
55 mat.template triangularView<UploA>() =
A;
58 mat -= (
B *
sigma).
template triangularView<UploA>();
60 mat -= (
B *
sigma).
template triangularView<UploB>().transpose();
62 fac.compute(
mat, UploA);
69 template <
int UploA,
int UploB>
70 class SymShiftInvertHelper<true, false, UploA, UploB>
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)
76 using Matrix =
typename ArgB::PlainObject;
79 mat.template triangularView<UploB>() = -
sigma *
B;
82 mat +=
A.template triangularView<UploB>();
84 mat +=
A.template triangularView<UploA>().transpose();
86 fac.compute(
mat, UploB);
126 typename StorageIndexA =
int,
typename StorageIndexB =
int>
142 using ASparse = std::is_same<TypeA, Eigen::Sparse>;
150 using BSparse = std::is_same<TypeB, Eigen::Sparse>;
167 using FacType =
typename std::conditional<
190 template <
typename DerivedA,
typename DerivedB>
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)");
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)");
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");
222 using Helper = SymShiftInvertHelper<AIsSparse, BIsSparse, UploA, UploB>;
225 throw std::invalid_argument(
"SymShiftInvert: factorization failed with the given shift");
245 #endif // SPECTRA_SYM_SHIFT_INVERT_H
ConstGenericMatrixA m_matA
@ Successful
Computation was successful.
typename std::conditional< ASparse::value, MatrixB, MatrixA >::type DenseType
typename std::conditional< ASparse::value, SparseTypeA, DenseTypeA >::type MatrixA
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
void perform_op(const Scalar *x_in, Scalar *y_out) const
static const double sigma
typename std::conditional< ASparse::value &&BSparse::value, MatrixA, DenseType >::type ResType
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > Matrix
void set_shift(const Scalar &sigma)
std::is_same< Eigen::Sparse, Eigen::Sparse > BSparse
A matrix or vector expression mapping an existing array of data.
typename std::conditional< ASparse::value &&BSparse::value, Eigen::SparseLU< ResType >, BKLDLT< Scalar > >::type FacType
Eigen::SparseMatrix< double > Sparse
ConstGenericMatrixB m_matB
SymShiftInvert(const Eigen::EigenBase< DerivedA > &A, const Eigen::EigenBase< DerivedB > &B)
The matrix class, also used for vectors and row-vectors.
Sparse supernodal LU factorization for general matrices.
std::is_same< Eigen::Sparse, Eigen::Sparse > ASparse
typename std::conditional< BSparse::value, SparseTypeB, DenseTypeB >::type MatrixB
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
gtsam
Author(s):
autogenerated on Thu Apr 10 2025 03:04:25