mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-06-04 18:54:00 +08:00
Add condition estimation to Cholesky (LLT) factorization.
This commit is contained in:
parent
fb8dccc23e
commit
f54137606e
@ -135,6 +135,16 @@ template<typename _MatrixType, int _UpLo> class LLT
|
|||||||
template<typename InputType>
|
template<typename InputType>
|
||||||
LLT& compute(const EigenBase<InputType>& matrix);
|
LLT& compute(const EigenBase<InputType>& matrix);
|
||||||
|
|
||||||
|
/** \returns an estimate of the reciprocal condition number of the matrix of
|
||||||
|
* which *this is the Cholesky decomposition.
|
||||||
|
*/
|
||||||
|
RealScalar rcond() const
|
||||||
|
{
|
||||||
|
eigen_assert(m_isInitialized && "LLT is not initialized.");
|
||||||
|
eigen_assert(m_info == Success && "LLT failed because matrix appears to be negative");
|
||||||
|
return ConditionEstimator<LLT<MatrixType, UpLo>, true >::rcond(m_l1_norm, *this);
|
||||||
|
}
|
||||||
|
|
||||||
/** \returns the LLT decomposition matrix
|
/** \returns the LLT decomposition matrix
|
||||||
*
|
*
|
||||||
* TODO: document the storage layout
|
* TODO: document the storage layout
|
||||||
@ -183,6 +193,7 @@ template<typename _MatrixType, int _UpLo> class LLT
|
|||||||
* The strict upper part is not used and even not initialized.
|
* The strict upper part is not used and even not initialized.
|
||||||
*/
|
*/
|
||||||
MatrixType m_matrix;
|
MatrixType m_matrix;
|
||||||
|
RealScalar m_l1_norm;
|
||||||
bool m_isInitialized;
|
bool m_isInitialized;
|
||||||
ComputationInfo m_info;
|
ComputationInfo m_info;
|
||||||
};
|
};
|
||||||
@ -393,6 +404,26 @@ LLT<MatrixType,_UpLo>& LLT<MatrixType,_UpLo>::compute(const EigenBase<InputType>
|
|||||||
m_matrix.resize(size, size);
|
m_matrix.resize(size, size);
|
||||||
m_matrix = a.derived();
|
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_isInitialized = true;
|
m_isInitialized = true;
|
||||||
bool ok = Traits::inplace_decomposition(m_matrix);
|
bool ok = Traits::inplace_decomposition(m_matrix);
|
||||||
m_info = ok ? Success : NumericalIssue;
|
m_info = ok ? Success : NumericalIssue;
|
||||||
|
@ -13,11 +13,11 @@
|
|||||||
namespace Eigen {
|
namespace Eigen {
|
||||||
|
|
||||||
namespace internal {
|
namespace internal {
|
||||||
template <typename Decomposition, bool IsComplex>
|
template <typename Decomposition, bool IsSelfAdjoint, bool IsComplex>
|
||||||
struct EstimateInverseMatrixL1NormImpl {};
|
struct EstimateInverseMatrixL1NormImpl {};
|
||||||
} // namespace internal
|
} // namespace internal
|
||||||
|
|
||||||
template <typename Decomposition>
|
template <typename Decomposition, bool IsSelfAdjoint = false>
|
||||||
class ConditionEstimator {
|
class ConditionEstimator {
|
||||||
public:
|
public:
|
||||||
typedef typename Decomposition::MatrixType MatrixType;
|
typedef typename Decomposition::MatrixType MatrixType;
|
||||||
@ -101,7 +101,8 @@ class ConditionEstimator {
|
|||||||
return 0;
|
return 0;
|
||||||
}
|
}
|
||||||
return internal::EstimateInverseMatrixL1NormImpl<
|
return internal::EstimateInverseMatrixL1NormImpl<
|
||||||
Decomposition, NumTraits<Scalar>::IsComplex>::compute(dec);
|
Decomposition, IsSelfAdjoint,
|
||||||
|
NumTraits<Scalar>::IsComplex != 0>::compute(dec);
|
||||||
}
|
}
|
||||||
|
|
||||||
/**
|
/**
|
||||||
@ -116,9 +117,27 @@ class ConditionEstimator {
|
|||||||
|
|
||||||
namespace internal {
|
namespace internal {
|
||||||
|
|
||||||
|
template <typename Decomposition, typename Vector, bool IsSelfAdjoint = false>
|
||||||
|
struct solve_helper {
|
||||||
|
static inline Vector solve_adjoint(const Decomposition& dec,
|
||||||
|
const Vector& v) {
|
||||||
|
return dec.adjoint().solve(v);
|
||||||
|
}
|
||||||
|
};
|
||||||
|
|
||||||
|
// Partial specialization for self_adjoint matrices.
|
||||||
|
template <typename Decomposition, typename Vector>
|
||||||
|
struct solve_helper<Decomposition, Vector, true> {
|
||||||
|
static inline Vector solve_adjoint(const Decomposition& dec,
|
||||||
|
const Vector& v) {
|
||||||
|
return dec.solve(v);
|
||||||
|
}
|
||||||
|
};
|
||||||
|
|
||||||
|
|
||||||
// Partial specialization for real matrices.
|
// Partial specialization for real matrices.
|
||||||
template <typename Decomposition>
|
template <typename Decomposition, bool IsSelfAdjoint>
|
||||||
struct EstimateInverseMatrixL1NormImpl<Decomposition, 0> {
|
struct EstimateInverseMatrixL1NormImpl<Decomposition, IsSelfAdjoint, false> {
|
||||||
typedef typename Decomposition::MatrixType MatrixType;
|
typedef typename Decomposition::MatrixType MatrixType;
|
||||||
typedef typename internal::traits<MatrixType>::Scalar Scalar;
|
typedef typename internal::traits<MatrixType>::Scalar Scalar;
|
||||||
typedef typename internal::plain_col_type<MatrixType>::type Vector;
|
typedef typename internal::plain_col_type<MatrixType>::type Vector;
|
||||||
@ -152,7 +171,7 @@ struct EstimateInverseMatrixL1NormImpl<Decomposition, 0> {
|
|||||||
int old_v_max_abs_index = v_max_abs_index;
|
int old_v_max_abs_index = v_max_abs_index;
|
||||||
for (int k = 0; k < 4; ++k) {
|
for (int k = 0; k < 4; ++k) {
|
||||||
// argmax |inv(matrix)^T * sign_vector|
|
// argmax |inv(matrix)^T * sign_vector|
|
||||||
v = dec.adjoint().solve(sign_vector);
|
v = solve_helper<Decomposition, Vector, IsSelfAdjoint>::solve_adjoint(dec, sign_vector);
|
||||||
v.cwiseAbs().maxCoeff(&v_max_abs_index);
|
v.cwiseAbs().maxCoeff(&v_max_abs_index);
|
||||||
if (v_max_abs_index == old_v_max_abs_index) {
|
if (v_max_abs_index == old_v_max_abs_index) {
|
||||||
// Break if the solution stagnated.
|
// Break if the solution stagnated.
|
||||||
@ -200,8 +219,8 @@ struct EstimateInverseMatrixL1NormImpl<Decomposition, 0> {
|
|||||||
};
|
};
|
||||||
|
|
||||||
// Partial specialization for complex matrices.
|
// Partial specialization for complex matrices.
|
||||||
template <typename Decomposition>
|
template <typename Decomposition, bool IsSelfAdjoint>
|
||||||
struct EstimateInverseMatrixL1NormImpl<Decomposition, 1> {
|
struct EstimateInverseMatrixL1NormImpl<Decomposition, IsSelfAdjoint, true> {
|
||||||
typedef typename Decomposition::MatrixType MatrixType;
|
typedef typename Decomposition::MatrixType MatrixType;
|
||||||
typedef typename internal::traits<MatrixType>::Scalar Scalar;
|
typedef typename internal::traits<MatrixType>::Scalar Scalar;
|
||||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||||
@ -238,7 +257,7 @@ struct EstimateInverseMatrixL1NormImpl<Decomposition, 1> {
|
|||||||
RealVector abs_v = v.cwiseAbs();
|
RealVector abs_v = v.cwiseAbs();
|
||||||
const Vector psi =
|
const Vector psi =
|
||||||
(abs_v.array() == 0).select(v.cwiseQuotient(abs_v), ones);
|
(abs_v.array() == 0).select(v.cwiseQuotient(abs_v), ones);
|
||||||
v = dec.adjoint().solve(psi);
|
v = solve_helper<Decomposition, Vector, IsSelfAdjoint>::solve_adjoint(dec, psi);
|
||||||
const RealVector z = v.real();
|
const RealVector z = v.real();
|
||||||
z.cwiseAbs().maxCoeff(&v_max_abs_index);
|
z.cwiseAbs().maxCoeff(&v_max_abs_index);
|
||||||
if (v_max_abs_index == old_v_max_abs_index) {
|
if (v_max_abs_index == old_v_max_abs_index) {
|
||||||
|
@ -17,6 +17,12 @@
|
|||||||
#include <Eigen/Cholesky>
|
#include <Eigen/Cholesky>
|
||||||
#include <Eigen/QR>
|
#include <Eigen/QR>
|
||||||
|
|
||||||
|
template<typename MatrixType, int UpLo>
|
||||||
|
typename MatrixType::RealScalar matrix_l1_norm(const MatrixType& m) {
|
||||||
|
MatrixType symm = m.template selfadjointView<UpLo>();
|
||||||
|
return symm.cwiseAbs().colwise().sum().maxCoeff();
|
||||||
|
}
|
||||||
|
|
||||||
template<typename MatrixType,template <typename,int> class CholType> void test_chol_update(const MatrixType& symm)
|
template<typename MatrixType,template <typename,int> class CholType> void test_chol_update(const MatrixType& symm)
|
||||||
{
|
{
|
||||||
typedef typename MatrixType::Scalar Scalar;
|
typedef typename MatrixType::Scalar Scalar;
|
||||||
@ -85,6 +91,14 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
|
|||||||
matX = chollo.solve(matB);
|
matX = chollo.solve(matB);
|
||||||
VERIFY_IS_APPROX(symm * matX, 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 = chollo.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 = chollo.rcond();
|
||||||
|
VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10);
|
||||||
|
|
||||||
// test the upper mode
|
// test the upper mode
|
||||||
LLT<SquareMatrixType,Upper> cholup(symmUp);
|
LLT<SquareMatrixType,Upper> cholup(symmUp);
|
||||||
VERIFY_IS_APPROX(symm, cholup.reconstructedMatrix());
|
VERIFY_IS_APPROX(symm, cholup.reconstructedMatrix());
|
||||||
@ -93,6 +107,15 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
|
|||||||
matX = cholup.solve(matB);
|
matX = cholup.solve(matB);
|
||||||
VERIFY_IS_APPROX(symm * matX, 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 = cholup.solve(MatrixType::Identity(rows,cols));
|
||||||
|
rcond = (RealScalar(1) / matrix_l1_norm<MatrixType, Upper>(symmUp)) /
|
||||||
|
matrix_l1_norm<MatrixType, Upper>(symmUp_inverse);
|
||||||
|
rcond_est = cholup.rcond();
|
||||||
|
VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10);
|
||||||
|
|
||||||
|
|
||||||
MatrixType neg = -symmLo;
|
MatrixType neg = -symmLo;
|
||||||
chollo.compute(neg);
|
chollo.compute(neg);
|
||||||
VERIFY(chollo.info()==NumericalIssue);
|
VERIFY(chollo.info()==NumericalIssue);
|
||||||
|
10
test/lu.cpp
10
test/lu.cpp
@ -151,10 +151,11 @@ template<typename MatrixType> void lu_invertible()
|
|||||||
MatrixType m1_inverse = lu.inverse();
|
MatrixType m1_inverse = lu.inverse();
|
||||||
VERIFY_IS_APPROX(m2, m1_inverse*m3);
|
VERIFY_IS_APPROX(m2, m1_inverse*m3);
|
||||||
|
|
||||||
// Test condition number estimation.
|
// Verify that the estimated condition number is within a factor of 10 of the
|
||||||
|
// truth.
|
||||||
RealScalar rcond = (RealScalar(1) / matrix_l1_norm(m1)) / matrix_l1_norm(m1_inverse);
|
RealScalar rcond = (RealScalar(1) / matrix_l1_norm(m1)) / matrix_l1_norm(m1_inverse);
|
||||||
// Verify that the estimate is within a factor of 10 of the truth.
|
const RealScalar rcond_est = lu.rcond();
|
||||||
VERIFY(lu.rcond() > rcond / 10 && lu.rcond() < rcond * 10);
|
VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10);
|
||||||
|
|
||||||
// test solve with transposed
|
// test solve with transposed
|
||||||
lu.template _solve_impl_transposed<false>(m3, m2);
|
lu.template _solve_impl_transposed<false>(m3, m2);
|
||||||
@ -199,7 +200,8 @@ template<typename MatrixType> void lu_partial_piv()
|
|||||||
// Test condition number estimation.
|
// Test condition number estimation.
|
||||||
RealScalar rcond = (RealScalar(1) / matrix_l1_norm(m1)) / matrix_l1_norm(m1_inverse);
|
RealScalar rcond = (RealScalar(1) / matrix_l1_norm(m1)) / matrix_l1_norm(m1_inverse);
|
||||||
// Verify that the estimate is within a factor of 10 of the truth.
|
// Verify that the estimate is within a factor of 10 of the truth.
|
||||||
VERIFY(plu.rcond() > rcond / 10 && plu.rcond() < rcond * 10);
|
const RealScalar rcond_est = plu.rcond();
|
||||||
|
VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10);
|
||||||
|
|
||||||
// test solve with transposed
|
// test solve with transposed
|
||||||
plu.template _solve_impl_transposed<false>(m3, m2);
|
plu.template _solve_impl_transposed<false>(m3, m2);
|
||||||
|
Loading…
x
Reference in New Issue
Block a user