11 #ifndef EIGEN_TRIDIAGONALIZATION_H 12 #define EIGEN_TRIDIAGONALIZATION_H 18 template<
typename MatrixType>
struct TridiagonalizationMatrixTReturnType;
19 template<
typename MatrixType>
20 struct traits<TridiagonalizationMatrixTReturnType<MatrixType> >
21 :
public traits<typename MatrixType::PlainObject>
23 typedef typename MatrixType::PlainObject ReturnType;
27 template<
typename MatrixType,
typename CoeffVectorType>
28 void tridiagonalization_inplace(MatrixType& matA, CoeffVectorType& hCoeffs);
70 typedef typename MatrixType::Scalar Scalar;
75 Size = MatrixType::RowsAtCompileTime,
76 SizeMinusOne = Size ==
Dynamic ?
Dynamic : (Size > 1 ? Size - 1 : 1),
77 Options = MatrixType::Options,
78 MaxSize = MatrixType::MaxRowsAtCompileTime,
79 MaxSizeMinusOne = MaxSize ==
Dynamic ?
Dynamic : (MaxSize > 1 ? MaxSize - 1 : 1)
83 typedef typename internal::plain_col_type<MatrixType, RealScalar>::type DiagonalType;
85 typedef typename internal::remove_all<typename MatrixType::RealReturnType>::type MatrixTypeRealView;
86 typedef internal::TridiagonalizationMatrixTReturnType<MatrixTypeRealView> MatrixTReturnType;
88 typedef typename internal::conditional<NumTraits<Scalar>::IsComplex,
89 typename internal::add_const_on_value_type<typename Diagonal<const MatrixType>::RealReturnType>::type,
91 >::type DiagonalReturnType;
93 typedef typename internal::conditional<NumTraits<Scalar>::IsComplex,
94 typename internal::add_const_on_value_type<
typename Diagonal<
const MatrixType, -1>::RealReturnType>::type,
96 >::type SubDiagonalReturnType;
114 : m_matrix(size,size),
115 m_hCoeffs(size > 1 ? size-1 : 1),
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)
135 internal::tridiagonalization_inplace(m_matrix, m_hCoeffs);
136 m_isInitialized =
true;
156 template<
typename InputType>
160 m_hCoeffs.resize(matrix.
rows()-1, 1);
161 internal::tridiagonalization_inplace(m_matrix, m_hCoeffs);
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;
301 CoeffVectorType m_hCoeffs;
302 bool m_isInitialized;
305 template<
typename MatrixType>
306 typename Tridiagonalization<MatrixType>::DiagonalReturnType
309 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
310 return m_matrix.diagonal().real();
313 template<
typename MatrixType>
314 typename Tridiagonalization<MatrixType>::SubDiagonalReturnType
317 eigen_assert(m_isInitialized &&
"Tridiagonalization is not initialized.");
318 return m_matrix.template diagonal<-1>().
real();
346 template<
typename MatrixType,
typename CoeffVectorType>
350 typedef typename MatrixType::Scalar Scalar;
351 typedef typename MatrixType::RealScalar RealScalar;
352 Index n = matA.rows();
353 eigen_assert(n==matA.cols());
354 eigen_assert(n==hCoeffs.size()+1 || n==1);
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));
382 int Size=MatrixType::ColsAtCompileTime,
384 struct tridiagonalization_inplace_selector;
426 template<
typename MatrixType,
typename DiagonalType,
typename SubDiagonalType>
429 eigen_assert(mat.cols()==mat.rows() && diag.size()==mat.rows() && subdiag.
size()==mat.rows()-1);
430 tridiagonalization_inplace_selector<MatrixType>::run(mat, diag, subdiag, extractQ);
436 template<
typename MatrixType,
int Size,
bool IsComplex>
437 struct tridiagonalization_inplace_selector
441 template<
typename DiagonalType,
typename SubDiagonalType>
444 CoeffVectorType hCoeffs(mat.cols()-1);
445 tridiagonalization_inplace(mat,hCoeffs);
446 diag = mat.diagonal().real();
447 subdiag = mat.template diagonal<-1>().
real();
449 mat = HouseholderSequenceType(mat, hCoeffs.conjugate())
450 .setLength(mat.rows() - 1)
459 template<
typename MatrixType>
460 struct tridiagonalization_inplace_selector<MatrixType,3,false>
462 typedef typename MatrixType::Scalar Scalar;
463 typedef typename MatrixType::RealScalar RealScalar;
465 template<
typename DiagonalType,
typename SubDiagonalType>
466 static void run(MatrixType& mat, DiagonalType& diag,
SubDiagonalType& subdiag,
bool extractQ)
469 const RealScalar tol = (std::numeric_limits<RealScalar>::min)();
471 RealScalar v1norm2 = numext::abs2(mat(2,0));
476 subdiag[0] = mat(1,0);
477 subdiag[1] = mat(2,1);
483 RealScalar beta =
sqrt(numext::abs2(mat(1,0)) + v1norm2);
484 RealScalar invBeta = RealScalar(1)/beta;
485 Scalar m01 = mat(1,0) * invBeta;
486 Scalar m02 = mat(2,0) * invBeta;
487 Scalar q = RealScalar(2)*m01*mat(2,1) + m02*(mat(2,2) - mat(1,1));
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>
506 struct tridiagonalization_inplace_selector<MatrixType,1,IsComplex>
508 typedef typename MatrixType::Scalar Scalar;
510 template<
typename DiagonalType,
typename SubDiagonalType>
511 static void run(MatrixType& mat, DiagonalType& diag,
SubDiagonalType&,
bool extractQ)
513 diag(0,0) = numext::real(mat(0,0));
515 mat(0,0) = Scalar(1);
526 template<
typename MatrixType>
struct TridiagonalizationMatrixTReturnType
527 :
public ReturnByValue<TridiagonalizationMatrixTReturnType<MatrixType> >
534 TridiagonalizationMatrixTReturnType(
const MatrixType& mat) : m_matrix(mat) { }
536 template <
typename ResultType>
537 inline void evalTo(ResultType& result)
const 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(); }
549 typename MatrixType::Nested m_matrix;
556 #endif // EIGEN_TRIDIAGONALIZATION_H Tridiagonalization(Index size=Size==Dynamic ? 2 :Size)
Default constructor.
Definition: Tridiagonalization.h:113
DiagonalReturnType diagonal() const
Returns the diagonal of the tridiagonal matrix T in the decomposition.
Definition: Tridiagonalization.h:307
HouseholderSequence< MatrixType, typename internal::remove_all< typename CoeffVectorType::ConjugateReturnType >::type > HouseholderSequenceType
Return type of matrixQ()
Definition: Tridiagonalization.h:99
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Eigen::Index Index
Definition: Tridiagonalization.h:72
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
HouseholderSequenceType matrixQ() const
Returns the unitary matrix Q in the decomposition.
Definition: Tridiagonalization.h:240
Namespace containing all symbols from the Eigen library.
Definition: Core:287
Tridiagonalization & compute(const EigenBase< InputType > &matrix)
Computes tridiagonal decomposition of given matrix.
Definition: Tridiagonalization.h:157
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:150
Tridiagonal decomposition of a selfadjoint matrix.
Definition: Tridiagonalization.h:63
Derived & derived()
Definition: EigenBase.h:45
Definition: EigenBase.h:29
Sequence of Householder reflections acting on subspaces with decreasing size.
Definition: ForwardDeclarations.h:262
Scalar & coeffRef(Index rowId, Index colId)
Definition: PlainObjectBase.h:183
const MatrixType & packedMatrix() const
Returns the internal representation of the decomposition.
Definition: Tridiagonalization.h:219
Index size() const
Definition: EigenBase.h:66
CoeffVectorType householderCoefficients() const
Returns the Householder coefficients.
Definition: Tridiagonalization.h:182
Index rows() const
Definition: EigenBase.h:59
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
MatrixTReturnType matrixT() const
Returns an expression of the tridiagonal matrix T in the decomposition.
Definition: Tridiagonalization.h:265
Tridiagonalization(const EigenBase< InputType > &matrix)
Constructor; computes tridiagonal decomposition of given matrix.
Definition: Tridiagonalization.h:130
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_real_op< typename Derived::Scalar >, const Derived > real(const Eigen::ArrayBase< Derived > &x)
Definition: Eigen_Colamd.h:50
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:63
const int Dynamic
Definition: Constants.h:21
_MatrixType MatrixType
Synonym for the template parameter _MatrixType.
Definition: Tridiagonalization.h:68
SubDiagonalReturnType subDiagonal() const
Returns the subdiagonal of the tridiagonal matrix T in the decomposition.
Definition: Tridiagonalization.h:315