11 #ifndef EIGEN_TRIDIAGONALIZATION_H 12 #define EIGEN_TRIDIAGONALIZATION_H 19 template<
typename MatrixType>
21 :
public traits<typename MatrixType::PlainObject>
27 template<
typename MatrixType,
typename CoeffVectorType>
75 Size = MatrixType::RowsAtCompileTime,
77 Options = MatrixType::Options,
78 MaxSize = MatrixType::MaxRowsAtCompileTime,
79 MaxSizeMinusOne = MaxSize ==
Dynamic ?
Dynamic : (MaxSize > 1 ? MaxSize - 1 : 1)
93 typedef typename internal::conditional<NumTraits<Scalar>::IsComplex,
116 m_isInitialized(false)
129 template<
typename InputType>
131 : m_matrix(matrix.derived()),
132 m_hCoeffs(matrix.
cols() > 1 ? matrix.
cols()-1 : 1),
133 m_isInitialized(false)
136 m_isInitialized =
true;
156 template<
typename InputType>
160 m_hCoeffs.resize(matrix.
rows()-1, 1);
162 m_isInitialized =
true;
184 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
221 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
242 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
243 return HouseholderSequenceType(m_matrix, m_hCoeffs.conjugate())
244 .setLength(m_matrix.rows() - 1)
267 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
268 return MatrixTReturnType(m_matrix.real());
284 DiagonalReturnType
diagonal()
const;
296 SubDiagonalReturnType subDiagonal()
const;
305 template<
typename MatrixType>
309 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
310 return m_matrix.diagonal().real();
313 template<
typename MatrixType>
317 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
346 template<
typename MatrixType,
typename CoeffVectorType>
352 Index
n = matA.rows();
356 for (Index
i = 0;
i<n-1; ++
i)
358 Index remainingSize = n-
i-1;
361 matA.col(i).tail(remainingSize).makeHouseholderInPlace(h, beta);
365 matA.col(i).coeffRef(i+1) = 1;
367 hCoeffs.tail(n-i-1).noalias() = (matA.bottomRightCorner(remainingSize,remainingSize).template selfadjointView<Lower>()
368 * (
conj(h) * matA.col(i).tail(remainingSize)));
370 hCoeffs.tail(n-i-1) += (
conj(h)*
RealScalar(-0.5)*(hCoeffs.tail(remainingSize).dot(matA.col(i).tail(remainingSize)))) * matA.col(i).tail(n-i-1);
372 matA.bottomRightCorner(remainingSize, remainingSize).template selfadjointView<Lower>()
373 .rankUpdate(matA.col(i).tail(remainingSize), hCoeffs.tail(remainingSize),
Scalar(-1));
375 matA.col(i).coeffRef(i+1) = beta;
376 hCoeffs.coeffRef(i) =
h;
382 int Size=MatrixType::ColsAtCompileTime,
426 template<
typename MatrixType,
typename DiagonalType,
typename SubDiagonalType>
429 eigen_assert(mat.cols()==mat.rows() && diag.size()==mat.rows() && subdiag.size()==mat.rows()-1);
436 template<
typename MatrixType,
int Size,
bool IsComplex>
441 template<
typename DiagonalType,
typename SubDiagonalType>
442 static void run(MatrixType&
mat, DiagonalType&
diag, SubDiagonalType& subdiag,
bool extractQ)
444 CoeffVectorType hCoeffs(mat.cols()-1);
446 diag = mat.diagonal().real();
449 mat = HouseholderSequenceType(mat, hCoeffs.conjugate())
450 .setLength(mat.rows() - 1)
459 template<
typename MatrixType>
465 template<
typename DiagonalType,
typename SubDiagonalType>
466 static void run(MatrixType&
mat, DiagonalType&
diag, SubDiagonalType& subdiag,
bool extractQ)
476 subdiag[0] =
mat(1,0);
477 subdiag[1] =
mat(2,1);
485 Scalar m01 =
mat(1,0) * invBeta;
486 Scalar m02 =
mat(2,0) * invBeta;
488 diag[1] =
mat(1,1) + m02*
q;
489 diag[2] =
mat(2,2) - m02*
q;
491 subdiag[1] =
mat(2,1) - m01 *
q;
505 template<
typename MatrixType,
bool IsComplex>
510 template<
typename DiagonalType,
typename SubDiagonalType>
511 static void run(MatrixType&
mat, DiagonalType&
diag, SubDiagonalType&,
bool extractQ)
527 :
public ReturnByValue<TridiagonalizationMatrixTReturnType<MatrixType> >
536 template <
typename ResultType>
540 result.template diagonal<1>() = m_matrix.template diagonal<-1>().conjugate();
541 result.diagonal() = m_matrix.diagonal();
542 result.template
diagonal<-1>() = m_matrix.template diagonal<-1>();
545 Index
rows()
const {
return m_matrix.rows(); }
546 Index
cols()
const {
return m_matrix.cols(); }
556 #endif // EIGEN_TRIDIAGONALIZATION_H
HouseholderSequence< MatrixType, typename internal::remove_all< typename CoeffVectorType::ConjugateReturnType >::type > HouseholderSequenceType
Return type of matrixQ()
Matrix< Scalar, SizeMinusOne, 1, Options &~RowMajor, MaxSizeMinusOne, 1 > CoeffVectorType
const AutoDiffScalar< DerType > & conj(const AutoDiffScalar< DerType > &x)
Tridiagonalization< MatrixType >::HouseholderSequenceType HouseholderSequenceType
internal::remove_all< typename MatrixType::RealReturnType >::type MatrixTypeRealView
Matrix diag(const std::vector< Matrix > &Hs)
static void run(MatrixType &mat, DiagonalType &diag, SubDiagonalType &subdiag, bool extractQ)
HouseholderSequenceType matrixQ() const
Returns the unitary matrix Q in the decomposition.
internal::conditional< NumTraits< Scalar >::IsComplex, typename internal::add_const_on_value_type< typename Diagonal< const MatrixType,-1 >::RealReturnType >::type, const Diagonal< const MatrixType,-1 > >::type SubDiagonalReturnType
static void run(MatrixType &mat, DiagonalType &diag, SubDiagonalType &, bool extractQ)
void diagonal(const MatrixType &m)
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Namespace containing all symbols from the Eigen library.
Tridiagonalization & compute(const EigenBase< InputType > &matrix)
Computes tridiagonal decomposition of given matrix.
MatrixTReturnType matrixT() const
Returns an expression of the tridiagonal matrix T in the decomposition.
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Tridiagonal decomposition of a selfadjoint matrix.
MatrixType::Nested m_matrix
Sequence of Householder reflections acting on subspaces with decreasing size.
const MatrixType & packedMatrix() const
Returns the internal representation of the decomposition.
Tridiagonalization< MatrixType >::CoeffVectorType CoeffVectorType
CoeffVectorType m_hCoeffs
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
DiagonalReturnType diagonal() const
Returns the diagonal of the tridiagonal matrix T in the decomposition.
SubDiagonalReturnType subDiagonal() const
Returns the subdiagonal of the tridiagonal matrix T in the decomposition.
void evalTo(ResultType &result) const
EIGEN_DEVICE_FUNC const Scalar & q
MatrixType::Scalar Scalar
MatrixType::Scalar Scalar
Tridiagonalization(const EigenBase< InputType > &matrix)
Constructor; computes tridiagonal decomposition of given matrix.
NumTraits< Scalar >::Real RealScalar
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Abs2ReturnType abs2() const
void tridiagonalization_inplace(MatrixType &matA, CoeffVectorType &hCoeffs)
internal::TridiagonalizationMatrixTReturnType< MatrixTypeRealView > MatrixTReturnType
internal::plain_col_type< MatrixType, RealScalar >::type DiagonalType
TridiagonalizationMatrixTReturnType(const MatrixType &mat)
Constructor.
MatrixType::RealScalar RealScalar
CoeffVectorType householderCoefficients() const
Returns the Householder coefficients.
EIGEN_DEVICE_FUNC Index rows() const
internal::conditional< NumTraits< Scalar >::IsComplex, typename internal::add_const_on_value_type< typename Diagonal< const MatrixType >::RealReturnType >::type, const Diagonal< const MatrixType > >::type DiagonalReturnType
Tridiagonalization(Index size=Size==Dynamic?2:Size)
Default constructor.
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
MatrixType::Scalar Scalar
Matrix< RealScalar, SizeMinusOne, 1, Options &~RowMajor, MaxSizeMinusOne, 1 > SubDiagonalType
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
_MatrixType MatrixType
Synonym for the template parameter _MatrixType.
const AutoDiffScalar< DerType > & real(const AutoDiffScalar< DerType > &x)
NumTraits< Scalar >::Real RealScalar
EIGEN_DEVICE_FUNC Derived & derived()
MatrixType::PlainObject ReturnType
static void run(MatrixType &mat, DiagonalType &diag, SubDiagonalType &subdiag, bool extractQ)
ScalarWithExceptions conj(const ScalarWithExceptions &x)