11#ifndef EIGEN_TRIDIAGONALIZATION_H
12#define EIGEN_TRIDIAGONALIZATION_H
15#include "./InternalHeaderCheck.h"
21template <
typename MatrixType>
22struct TridiagonalizationMatrixTReturnType;
23template <
typename MatrixType>
24struct traits<TridiagonalizationMatrixTReturnType<MatrixType>> :
public traits<typename MatrixType::PlainObject> {
25 typedef typename MatrixType::PlainObject ReturnType;
29template <
typename MatrixType,
typename CoeffVectorType>
30EIGEN_DEVICE_FUNC
void tridiagonalization_inplace(MatrixType& matA, CoeffVectorType& hCoeffs);
65template <
typename MatrixType_>
71 typedef typename MatrixType::Scalar Scalar;
72 typedef typename NumTraits<Scalar>::Real RealScalar;
76 Size = MatrixType::RowsAtCompileTime,
77 SizeMinusOne = Size ==
Dynamic ?
Dynamic : (Size > 1 ? Size - 1 : 1),
78 Options = internal::traits<MatrixType>::Options,
79 MaxSize = MatrixType::MaxRowsAtCompileTime,
80 MaxSizeMinusOne = MaxSize ==
Dynamic ?
Dynamic : (MaxSize > 1 ? MaxSize - 1 : 1)
84 typedef typename internal::plain_col_type<MatrixType, RealScalar>::type DiagonalType;
86 typedef internal::remove_all_t<typename MatrixType::RealReturnType> MatrixTypeRealView;
87 typedef internal::TridiagonalizationMatrixTReturnType<MatrixTypeRealView> MatrixTReturnType;
89 typedef std::conditional_t<NumTraits<Scalar>::IsComplex,
90 internal::add_const_on_value_type_t<typename Diagonal<const MatrixType>::RealReturnType>,
94 typedef std::conditional_t<
96 internal::add_const_on_value_type_t<
typename Diagonal<
const MatrixType, -1>::RealReturnType>,
98 SubDiagonalReturnType;
117 : m_matrix(size, size), m_hCoeffs(size > 1 ? size - 1 : 1), m_isInitialized(false) {}
129 template <
typename InputType>
131 : m_matrix(matrix.derived()), m_hCoeffs(matrix.cols() > 1 ? matrix.cols() - 1 : 1), m_isInitialized(false) {
132 internal::tridiagonalization_inplace(m_matrix, m_hCoeffs);
133 m_isInitialized =
true;
153 template <
typename InputType>
157 internal::tridiagonalization_inplace(m_matrix, m_hCoeffs);
158 m_isInitialized =
true;
179 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
215 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
235 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
257 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
258 return MatrixTReturnType(m_matrix.real());
274 DiagonalReturnType
diagonal()
const;
290 CoeffVectorType m_hCoeffs;
291 bool m_isInitialized;
294template <
typename MatrixType>
296 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
297 return m_matrix.diagonal().real();
300template <
typename MatrixType>
302 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
303 return m_matrix.template diagonal<-1>().
real();
331template <
typename MatrixType,
typename CoeffVectorType>
332EIGEN_DEVICE_FUNC
void tridiagonalization_inplace(MatrixType& matA, CoeffVectorType& hCoeffs) {
334 typedef typename MatrixType::Scalar Scalar;
335 typedef typename MatrixType::RealScalar RealScalar;
336 Index n = matA.rows();
337 eigen_assert(n == matA.cols());
338 eigen_assert(n == hCoeffs.size() + 1 || n == 1);
340 for (
Index i = 0; i < n - 1; ++i) {
341 Index remainingSize = n - i - 1;
344 matA.col(i).tail(remainingSize).makeHouseholderInPlace(h, beta);
348 matA.col(i).coeffRef(i + 1) = 1;
350 hCoeffs.tail(n - i - 1).noalias() =
351 (matA.bottomRightCorner(remainingSize, remainingSize).template selfadjointView<Lower>() *
352 (
conj(h) * matA.col(i).tail(remainingSize)));
354 hCoeffs.tail(n - i - 1) +=
355 (
conj(h) * RealScalar(-0.5) * (hCoeffs.tail(remainingSize).dot(matA.col(i).tail(remainingSize)))) *
356 matA.col(i).tail(n - i - 1);
358 matA.bottomRightCorner(remainingSize, remainingSize)
359 .template selfadjointView<Lower>()
360 .rankUpdate(matA.col(i).tail(remainingSize), hCoeffs.tail(remainingSize), Scalar(-1));
362 matA.col(i).coeffRef(i + 1) = beta;
363 hCoeffs.coeffRef(i) = h;
368template <
typename MatrixType,
int Size = MatrixType::ColsAtCompileTime,
369 bool IsComplex = NumTraits<typename MatrixType::Scalar>::IsComplex>
370struct tridiagonalization_inplace_selector;
412template <
typename MatrixType,
typename DiagonalType,
typename SubDiagonalType,
typename CoeffVectorType,
413 typename WorkSpaceType>
414EIGEN_DEVICE_FUNC
void tridiagonalization_inplace(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag,
415 CoeffVectorType& hcoeffs, WorkSpaceType& workspace,
bool extractQ) {
416 eigen_assert(mat.cols() == mat.rows() && diag.size() == mat.rows() && subdiag.size() == mat.rows() - 1);
417 tridiagonalization_inplace_selector<MatrixType>::run(mat, diag, subdiag, hcoeffs, workspace, extractQ);
423template <
typename MatrixType,
int Size,
bool IsComplex>
424struct tridiagonalization_inplace_selector {
426 template <
typename DiagonalType,
typename SubDiagonalType,
typename CoeffVectorType,
typename WorkSpaceType>
427 static EIGEN_DEVICE_FUNC
void run(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag,
428 CoeffVectorType& hCoeffs, WorkSpaceType& workspace,
bool extractQ) {
429 tridiagonalization_inplace(mat, hCoeffs);
430 diag = mat.diagonal().real();
431 subdiag = mat.template diagonal<-1>().
real();
433 HouseholderSequenceType(mat, hCoeffs.conjugate()).setLength(mat.rows() - 1).setShift(1).evalTo(mat, workspace);
442template <
typename MatrixType>
443struct tridiagonalization_inplace_selector<MatrixType, 3, false> {
444 typedef typename MatrixType::Scalar Scalar;
445 typedef typename MatrixType::RealScalar RealScalar;
447 template <
typename DiagonalType,
typename SubDiagonalType,
typename CoeffVectorType,
typename WorkSpaceType>
448 static EIGEN_DEVICE_FUNC
void run(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag, CoeffVectorType&,
449 WorkSpaceType&,
bool extractQ) {
451 const RealScalar tol = (std::numeric_limits<RealScalar>::min)();
453 RealScalar v1norm2 = numext::abs2(mat(2, 0));
454 if (v1norm2 <= tol) {
457 subdiag[0] = mat(1, 0);
458 subdiag[1] = mat(2, 1);
459 if (extractQ) mat.setIdentity();
461 RealScalar beta =
sqrt(numext::abs2(mat(1, 0)) + v1norm2);
462 RealScalar invBeta = RealScalar(1) / beta;
463 Scalar m01 = mat(1, 0) * invBeta;
464 Scalar m02 = mat(2, 0) * invBeta;
465 Scalar q = RealScalar(2) * m01 * mat(2, 1) + m02 * (mat(2, 2) - mat(1, 1));
466 diag[1] = mat(1, 1) + m02 * q;
467 diag[2] = mat(2, 2) - m02 * q;
469 subdiag[1] = mat(2, 1) - m01 * q;
471 mat << 1, 0, 0, 0, m01, m02, 0, m02, -m01;
480template <
typename MatrixType,
bool IsComplex>
481struct tridiagonalization_inplace_selector<MatrixType, 1, IsComplex> {
482 typedef typename MatrixType::Scalar Scalar;
484 template <
typename DiagonalType,
typename SubDiagonalType,
typename CoeffVectorType,
typename WorkSpaceType>
485 static EIGEN_DEVICE_FUNC
void run(MatrixType& mat, DiagonalType& diag, SubDiagonalType&, CoeffVectorType&,
486 WorkSpaceType&,
bool extractQ) {
487 diag(0, 0) = numext::real(mat(0, 0));
488 if (extractQ) mat(0, 0) = Scalar(1);
499template <
typename MatrixType>
500struct TridiagonalizationMatrixTReturnType :
public ReturnByValue<TridiagonalizationMatrixTReturnType<MatrixType>> {
506 TridiagonalizationMatrixTReturnType(
const MatrixType& mat) : m_matrix(mat) {}
508 template <
typename ResultType>
509 inline void evalTo(ResultType& result)
const {
511 result.template diagonal<1>() = m_matrix.template diagonal<-1>().conjugate();
512 result.diagonal() = m_matrix.diagonal();
513 result.template diagonal<-1>() = m_matrix.template diagonal<-1>();
516 EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT {
return m_matrix.rows(); }
517 EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT {
return m_matrix.cols(); }
520 typename MatrixType::Nested m_matrix;
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition Diagonal.h:68
Sequence of Householder reflections acting on subspaces with decreasing size.
Definition HouseholderSequence.h:117
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
constexpr void resize(Index rows, Index cols)
Definition PlainObjectBase.h:294
Tridiagonal decomposition of a selfadjoint matrix.
Definition Tridiagonalization.h:66
Tridiagonalization(Index size=Size==Dynamic ? 2 :Size)
Default constructor.
Definition Tridiagonalization.h:116
DiagonalReturnType diagonal() const
Returns the diagonal of the tridiagonal matrix T in the decomposition.
Definition Tridiagonalization.h:295
MatrixType_ MatrixType
Synonym for the template parameter MatrixType_.
Definition Tridiagonalization.h:69
const MatrixType & packedMatrix() const
Returns the internal representation of the decomposition.
Definition Tridiagonalization.h:214
CoeffVectorType householderCoefficients() const
Returns the Householder coefficients.
Definition Tridiagonalization.h:178
HouseholderSequence< MatrixType, internal::remove_all_t< typename CoeffVectorType::ConjugateReturnType > > HouseholderSequenceType
Return type of matrixQ()
Definition Tridiagonalization.h:102
SubDiagonalReturnType subDiagonal() const
Returns the subdiagonal of the tridiagonal matrix T in the decomposition.
Definition Tridiagonalization.h:301
HouseholderSequenceType matrixQ() const
Returns the unitary matrix Q in the decomposition.
Definition Tridiagonalization.h:234
MatrixTReturnType matrixT() const
Returns an expression of the tridiagonal matrix T in the decomposition.
Definition Tridiagonalization.h:256
Tridiagonalization(const EigenBase< InputType > &matrix)
Constructor; computes tridiagonal decomposition of given matrix.
Definition Tridiagonalization.h:130
Tridiagonalization & compute(const EigenBase< InputType > &matrix)
Computes tridiagonal decomposition of given matrix.
Definition Tridiagonalization.h:154
Eigen::Index Index
Definition Tridiagonalization.h:73
Namespace containing all symbols from the Eigen library.
Definition Core:137
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_real_op< typename Derived::Scalar >, const Derived > real(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:83
const int Dynamic
Definition Constants.h:25
Definition EigenBase.h:33
Derived & derived()
Definition EigenBase.h:49
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition EigenBase.h:59
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition Meta.h:523