Spectra
1.2.0
Header-only C++ Library for Large Scale Eigenvalue Problems
Toggle main menu visibility
Loading...
Searching...
No Matches
DenseGenRealShiftSolve.h
1
// Copyright (C) 2016-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_DENSE_GEN_REAL_SHIFT_SOLVE_H
8
#define SPECTRA_DENSE_GEN_REAL_SHIFT_SOLVE_H
9
10
#include <Eigen/Core>
11
#include <Eigen/LU>
12
#include <stdexcept>
13
14
namespace
Spectra {
15
28
template
<
typename
Scalar_,
int
Flags = Eigen::ColMajor>
29
class
DenseGenRealShiftSolve
30
{
31
public
:
35
using
Scalar
= Scalar_;
36
37
private
:
38
using
Index = Eigen::Index;
39
using
Matrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Flags>;
40
using
Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
41
using
MapConstVec = Eigen::Map<const Vector>;
42
using
MapVec = Eigen::Map<Vector>;
43
using
ConstGenericMatrix =
const
Eigen::Ref<const Matrix>;
44
45
ConstGenericMatrix m_mat;
46
const
Index m_n;
47
Eigen::PartialPivLU<Matrix> m_solver;
48
49
public
:
58
template
<
typename
Derived>
59
DenseGenRealShiftSolve
(
const
Eigen::MatrixBase<Derived>& mat) :
60
m_mat(mat), m_n(mat.
rows
())
61
{
62
static_assert
(
63
static_cast<
int
>
(Derived::PlainObject::IsRowMajor) ==
static_cast<
int
>
(Matrix::IsRowMajor),
64
"DenseGenRealShiftSolve: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)"
);
65
66
if
(mat.rows() != mat.cols())
67
throw
std::invalid_argument(
"DenseGenRealShiftSolve: matrix must be square"
);
68
}
69
73
Index
rows
()
const
{
return
m_n; }
77
Index
cols
()
const
{
return
m_n; }
78
82
void
set_shift
(
const
Scalar
& sigma)
83
{
84
m_solver.compute(m_mat - sigma * Matrix::Identity(m_n, m_n));
85
}
86
93
// y_out = inv(A - sigma * I) * x_in
94
void
perform_op
(
const
Scalar
* x_in,
Scalar
* y_out)
const
95
{
96
MapConstVec x(x_in, m_n);
97
MapVec y(y_out, m_n);
98
y.noalias() = m_solver.solve(x);
99
}
100
};
101
102
}
// namespace Spectra
103
104
#endif
// SPECTRA_DENSE_GEN_REAL_SHIFT_SOLVE_H
Spectra::DenseGenRealShiftSolve::Scalar
Scalar_ Scalar
Definition
DenseGenRealShiftSolve.h:35
Spectra::DenseGenRealShiftSolve::rows
Index rows() const
Definition
DenseGenRealShiftSolve.h:73
Spectra::DenseGenRealShiftSolve::DenseGenRealShiftSolve
DenseGenRealShiftSolve(const Eigen::MatrixBase< Derived > &mat)
Definition
DenseGenRealShiftSolve.h:59
Spectra::DenseGenRealShiftSolve::set_shift
void set_shift(const Scalar &sigma)
Definition
DenseGenRealShiftSolve.h:82
Spectra::DenseGenRealShiftSolve::perform_op
void perform_op(const Scalar *x_in, Scalar *y_out) const
Definition
DenseGenRealShiftSolve.h:94
Spectra::DenseGenRealShiftSolve::cols
Index cols() const
Definition
DenseGenRealShiftSolve.h:77
Spectra
MatOp
DenseGenRealShiftSolve.h
Generated by
1.17.0