11 #ifndef UPPER_HESSENBERG_EIGEN_H 12 #define UPPER_HESSENBERG_EIGEN_H 15 #include <Eigen/Eigenvalues> 20 template <
typename Scalar =
double>
52 for (Index
j = 0;
j <
size; ++
j)
61 for (Index
n = size - 1;
n >= 0;
n--)
73 for (Index
i =
n - 1;
i >= 0;
i--)
76 Scalar r = m_matT.row(
i).segment(l,
n - l + 1).dot(m_matT.col(
n).segment(l,
n - l + 1));
98 Scalar t = (x * lastr - lastw * r) / denom;
103 m_matT.
coeffRef(
i + 1,
n) = (-lastr - y *
t) / lastw;
108 if ((eps * t) * t >
Scalar(1))
109 m_matT.col(
n).tail(size -
i) /=
t;
113 else if (q <
Scalar(0) &&
n > 0)
115 Scalar lastra(0), lastsa(0), lastw(0);
132 for (Index
i =
n - 2;
i >= 0;
i--)
134 Scalar ra = m_matT.row(
i).segment(l,
n - l + 1).dot(m_matT.col(
n - 1).segment(l,
n - l + 1));
135 Scalar sa = m_matT.row(
i).segment(l,
n - l + 1).dot(m_matT.col(
n).segment(l,
n - l + 1));
163 Complex cc =
Complex(x * lastra - lastw * ra + q * sa, x * lastsa - lastw * sa - q * ra) /
Complex(vr, vi);
181 if ((eps * t) * t >
Scalar(1))
182 m_matT.block(
i,
n - 1, size -
i, 2) /=
t;
193 for (Index
j = size - 1;
j >= 0;
j--)
195 m_tmp.noalias() = m_eivec.leftCols(
j + 1) * m_matT.col(
j).segment(0,
j + 1);
196 m_eivec.col(
j) = m_tmp;
202 m_n(0), m_computed(false)
206 m_n(mat.
rows()), m_computed(false)
216 if (mat.rows() != mat.cols())
217 throw std::invalid_argument(
"UpperHessenbergEigen: matrix must be square");
224 Matrix
Q = Matrix::Identity(m_n, m_n);
227 throw std::runtime_error(
"UpperHessenbergEigen: eigen decomposition failed");
229 m_matT = m_realSchur.
matrixT();
230 m_eivec = m_realSchur.
matrixU();
238 if (i == m_n - 1 || m_matT.
coeff(i + 1, i) ==
Scalar(0))
256 z = maxval *
sqrt(
abs(p0 * p0 + t0 * t1));
276 throw std::logic_error(
"UpperHessenbergEigen: need to call compute() first");
286 throw std::logic_error(
"UpperHessenbergEigen: need to call compute() first");
288 Index
n = m_eivec.
cols();
289 ComplexMatrix matV(n, n);
290 for (Index
j = 0;
j <
n; ++
j)
296 matV.col(
j) = m_eivec.col(
j).template cast<Complex>();
297 matV.col(
j).normalize();
302 for (Index
i = 0;
i <
n; ++
i)
307 matV.col(
j).normalize();
308 matV.col(
j + 1).normalize();
319 #endif // UPPER_HESSENBERG_EIGEN_H
Performs a real Schur decomposition of a square matrix.
void compute(ConstGenericMatrix &mat)
RealSchur & computeFromHessenberg(const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU)
Computes Schur decomposition of a Hessenberg matrix H = Z T Z^T.
const Eigen::Ref< const Matrix > ConstGenericMatrix
Eigen::Ref< Matrix > GenericMatrix
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
Eigen::Matrix< Complex, Eigen::Dynamic, Eigen::Dynamic > ComplexMatrix
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
static const Line3 l(Rot3(), 1, 1)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
std::complex< Scalar > Complex
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
UpperHessenbergEigen(ConstGenericMatrix &mat)
EIGEN_DEVICE_FUNC const Scalar & q
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
A matrix or vector expression mapping an existing expression.
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > Matrix
Eigen::Matrix< Complex, Eigen::Dynamic, 1 > ComplexVector
Eigen::RealSchur< Matrix > m_realSchur
const MatrixType & matrixU() const
Returns the orthogonal matrix in the Schur decomposition.
ComplexMatrix eigenvectors()
The quaternion class used to represent 3D orientations and rotations.
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 y set format x g set format y g set format x2 g set format y2 g set format z g set angles radians set nogrid set key title set key left top Right noreverse box linetype linewidth samplen spacing width set nolabel set noarrow set nologscale set logscale x set set pointsize set encoding default set nopolar set noparametric set set set set surface set nocontour set clabel set mapping cartesian set nohidden3d set cntrparam order set cntrparam linear set cntrparam levels auto set cntrparam points set size set set xzeroaxis lt lw set x2zeroaxis lt lw set yzeroaxis lt lw set y2zeroaxis lt lw set tics in set ticslevel set tics scale
EIGEN_DEVICE_FUNC const ImagReturnType imag() const
void doComputeEigenvectors()
Jet< T, N > sqrt(const Jet< T, N > &f)
ComputationInfo info() const
Reports whether previous computation was successful.
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
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
const ComplexVector & eigenvalues() const