mirror of
https://gitlab.com/libeigen/eigen.git
synced 2024-12-15 07:10:37 +08:00
Add rcond method to LDLT.
This commit is contained in:
parent
f54137606e
commit
9d51f7c457
@ -192,6 +192,15 @@ template<typename _MatrixType, int _UpLo> class LDLT
|
||||
template<typename InputType>
|
||||
LDLT& compute(const EigenBase<InputType>& matrix);
|
||||
|
||||
/** \returns an estimate of the reciprocal condition number of the matrix of
|
||||
* which *this is the LDLT decomposition.
|
||||
*/
|
||||
RealScalar rcond() const
|
||||
{
|
||||
eigen_assert(m_isInitialized && "LDLT is not initialized.");
|
||||
return ConditionEstimator<LDLT<MatrixType, UpLo>, true >::rcond(m_l1_norm, *this);
|
||||
}
|
||||
|
||||
template <typename Derived>
|
||||
LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1);
|
||||
|
||||
@ -241,6 +250,7 @@ template<typename _MatrixType, int _UpLo> class LDLT
|
||||
* is not stored), and the diagonal entries correspond to D.
|
||||
*/
|
||||
MatrixType m_matrix;
|
||||
RealScalar m_l1_norm;
|
||||
TranspositionType m_transpositions;
|
||||
TmpMatrixType m_temporary;
|
||||
internal::SignMatrix m_sign;
|
||||
@ -439,6 +449,26 @@ LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::compute(const EigenBase<InputTyp
|
||||
|
||||
m_matrix = a.derived();
|
||||
|
||||
// Compute matrix L1 norm = max abs column sum.
|
||||
m_l1_norm = RealScalar(0);
|
||||
if (_UpLo == Lower) {
|
||||
for (int col = 0; col < size; ++col) {
|
||||
const RealScalar abs_col_sum = m_matrix.col(col).tail(size - col).cwiseAbs().sum() +
|
||||
m_matrix.row(col).tail(col).cwiseAbs().sum();
|
||||
if (abs_col_sum > m_l1_norm) {
|
||||
m_l1_norm = abs_col_sum;
|
||||
}
|
||||
}
|
||||
} else {
|
||||
for (int col = 0; col < a.cols(); ++col) {
|
||||
const RealScalar abs_col_sum = m_matrix.col(col).tail(col).cwiseAbs().sum() +
|
||||
m_matrix.row(col).tail(size - col).cwiseAbs().sum();
|
||||
if (abs_col_sum > m_l1_norm) {
|
||||
m_l1_norm = abs_col_sum;
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
m_transpositions.resize(size);
|
||||
m_isInitialized = false;
|
||||
m_temporary.resize(size);
|
||||
|
@ -160,6 +160,15 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
|
||||
matX = ldltlo.solve(matB);
|
||||
VERIFY_IS_APPROX(symm * matX, matB);
|
||||
|
||||
// Verify that the estimated condition number is within a factor of 10 of the
|
||||
// truth.
|
||||
const MatrixType symmLo_inverse = ldltlo.solve(MatrixType::Identity(rows,cols));
|
||||
RealScalar rcond = (RealScalar(1) / matrix_l1_norm<MatrixType, Lower>(symmLo)) /
|
||||
matrix_l1_norm<MatrixType, Lower>(symmLo_inverse);
|
||||
RealScalar rcond_est = ldltlo.rcond();
|
||||
VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10);
|
||||
|
||||
|
||||
LDLT<SquareMatrixType,Upper> ldltup(symmUp);
|
||||
VERIFY_IS_APPROX(symm, ldltup.reconstructedMatrix());
|
||||
vecX = ldltup.solve(vecB);
|
||||
@ -167,6 +176,14 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
|
||||
matX = ldltup.solve(matB);
|
||||
VERIFY_IS_APPROX(symm * matX, matB);
|
||||
|
||||
// Verify that the estimated condition number is within a factor of 10 of the
|
||||
// truth.
|
||||
const MatrixType symmUp_inverse = ldltup.solve(MatrixType::Identity(rows,cols));
|
||||
rcond = (RealScalar(1) / matrix_l1_norm<MatrixType, Upper>(symmUp)) /
|
||||
matrix_l1_norm<MatrixType, Upper>(symmUp_inverse);
|
||||
rcond_est = ldltup.rcond();
|
||||
VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10);
|
||||
|
||||
VERIFY_IS_APPROX(MatrixType(ldltlo.matrixL().transpose().conjugate()), MatrixType(ldltlo.matrixU()));
|
||||
VERIFY_IS_APPROX(MatrixType(ldltlo.matrixU().transpose().conjugate()), MatrixType(ldltlo.matrixL()));
|
||||
VERIFY_IS_APPROX(MatrixType(ldltup.matrixL().transpose().conjugate()), MatrixType(ldltup.matrixU()));
|
||||
|
Loading…
Reference in New Issue
Block a user