aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/cholesky.cpp
diff options
context:
space:
mode:
authorGravatar Rasmus Munk Larsen <rmlarsen@google.com>2016-04-01 16:19:45 -0700
committerGravatar Rasmus Munk Larsen <rmlarsen@google.com>2016-04-01 16:19:45 -0700
commitf54137606eb6d68cbafd10d90013e254b26137ed (patch)
treebfb7752ec5bd7d6334c7fb139a898e148b203327 /test/cholesky.cpp
parentfb8dccc23e5f717319c230c2701a5fbf1d3c3975 (diff)
Add condition estimation to Cholesky (LLT) factorization.
Diffstat (limited to 'test/cholesky.cpp')
-rw-r--r--test/cholesky.cpp47
1 files changed, 35 insertions, 12 deletions
diff --git a/test/cholesky.cpp b/test/cholesky.cpp
index d652af5bf..8a21cdbd5 100644
--- a/test/cholesky.cpp
+++ b/test/cholesky.cpp
@@ -17,6 +17,12 @@
#include <Eigen/Cholesky>
#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)
{
typedef typename MatrixType::Scalar Scalar;
@@ -77,7 +83,7 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
{
SquareMatrixType symmUp = symm.template triangularView<Upper>();
SquareMatrixType symmLo = symm.template triangularView<Lower>();
-
+
LLT<SquareMatrixType,Lower> chollo(symmLo);
VERIFY_IS_APPROX(symm, chollo.reconstructedMatrix());
vecX = chollo.solve(vecB);
@@ -85,6 +91,14 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
matX = chollo.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 = 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
LLT<SquareMatrixType,Upper> cholup(symmUp);
VERIFY_IS_APPROX(symm, cholup.reconstructedMatrix());
@@ -93,6 +107,15 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
matX = cholup.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 = 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;
chollo.compute(neg);
VERIFY(chollo.info()==NumericalIssue);
@@ -101,7 +124,7 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
VERIFY_IS_APPROX(MatrixType(chollo.matrixU().transpose().conjugate()), MatrixType(chollo.matrixL()));
VERIFY_IS_APPROX(MatrixType(cholup.matrixL().transpose().conjugate()), MatrixType(cholup.matrixU()));
VERIFY_IS_APPROX(MatrixType(cholup.matrixU().transpose().conjugate()), MatrixType(cholup.matrixL()));
-
+
// test some special use cases of SelfCwiseBinaryOp:
MatrixType m1 = MatrixType::Random(rows,cols), m2(rows,cols);
m2 = m1;
@@ -167,7 +190,7 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
// restore
if(sign == -1)
symm = -symm;
-
+
// check matrices coming from linear constraints with Lagrange multipliers
if(rows>=3)
{
@@ -183,7 +206,7 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
vecX = ldltlo.solve(vecB);
VERIFY_IS_APPROX(A * vecX, vecB);
}
-
+
// check non-full rank matrices
if(rows>=3)
{
@@ -199,7 +222,7 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
vecX = ldltlo.solve(vecB);
VERIFY_IS_APPROX(A * vecX, vecB);
}
-
+
// check matrices with a wide spectrum
if(rows>=3)
{
@@ -225,7 +248,7 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
{
RealScalar large_tol = std::sqrt(test_precision<RealScalar>());
VERIFY((A * vecX).isApprox(vecB, large_tol));
-
+
++g_test_level;
VERIFY_IS_APPROX(A * vecX,vecB);
--g_test_level;
@@ -314,14 +337,14 @@ template<typename MatrixType> void cholesky_bug241(const MatrixType& m)
}
// LDLT is not guaranteed to work for indefinite matrices, but happens to work fine if matrix is diagonal.
-// This test checks that LDLT reports correctly that matrix is indefinite.
+// This test checks that LDLT reports correctly that matrix is indefinite.
// See http://forum.kde.org/viewtopic.php?f=74&t=106942 and bug 736
template<typename MatrixType> void cholesky_definiteness(const MatrixType& m)
{
eigen_assert(m.rows() == 2 && m.cols() == 2);
MatrixType mat;
LDLT<MatrixType> ldlt(2);
-
+
{
mat << 1, 0, 0, -1;
ldlt.compute(mat);
@@ -384,11 +407,11 @@ void test_cholesky()
CALL_SUBTEST_3( cholesky_definiteness(Matrix2d()) );
CALL_SUBTEST_4( cholesky(Matrix3f()) );
CALL_SUBTEST_5( cholesky(Matrix4d()) );
-
- s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE);
+
+ s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE);
CALL_SUBTEST_2( cholesky(MatrixXd(s,s)) );
TEST_SET_BUT_UNUSED_VARIABLE(s)
-
+
s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2);
CALL_SUBTEST_6( cholesky_cplx(MatrixXcd(s,s)) );
TEST_SET_BUT_UNUSED_VARIABLE(s)
@@ -402,6 +425,6 @@ void test_cholesky()
// Test problem size constructors
CALL_SUBTEST_9( LLT<MatrixXf>(10) );
CALL_SUBTEST_9( LDLT<MatrixXf>(10) );
-
+
TEST_SET_BUT_UNUSED_VARIABLE(nb_temporaries)
}