eigen/test/cholesky.cpp
Benoit Jacob a1ba995f05 Many improvements in LLT and LDLT:
* in LDLT, support the negative semidefinite case
* fix bad floating-point comparisons, improves greatly the accuracy of methods like
  isPositiveDefinite() and rank()
* simplifications
* identify (but not resolve) bug: claim that only triangular part is used, is inaccurate
* expanded unit-tests
2009-03-30 21:45:45 +00:00

200 lines
6.4 KiB
C++

// This file is part of Eigen, a lightweight C++ template library
// for linear algebra. Eigen itself is part of the KDE project.
//
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 3 of the License, or (at your option) any later version.
//
// Alternatively, you can redistribute it and/or
// modify it under the terms of the GNU General Public License as
// published by the Free Software Foundation; either version 2 of
// the License, or (at your option) any later version.
//
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.
#define EIGEN_NO_ASSERTION_CHECKING
#include "main.h"
#include <Eigen/Cholesky>
#include <Eigen/QR>
#ifdef HAS_GSL
#include "gsl_helper.h"
#endif
template<typename MatrixType> void cholesky(const MatrixType& m)
{
/* this test covers the following files:
LLT.h LDLT.h
*/
int rows = m.rows();
int cols = m.cols();
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
MatrixType a0 = MatrixType::Random(rows,cols);
VectorType vecB = VectorType::Random(rows), vecX(rows);
MatrixType matB = MatrixType::Random(rows,cols), matX(rows,cols);
SquareMatrixType symm = a0 * a0.adjoint();
// let's make sure the matrix is not singular or near singular
MatrixType a1 = MatrixType::Random(rows,cols);
symm += a1 * a1.adjoint();
// to test if really Cholesky only uses the upper triangular part, uncomment the following
// FIXME: currently that fails !!
//symm.template part<StrictlyLowerTriangular>().setZero();
#ifdef HAS_GSL
if (ei_is_same_type<RealScalar,double>::ret)
{
typedef GslTraits<Scalar> Gsl;
typename Gsl::Matrix gMatA=0, gSymm=0;
typename Gsl::Vector gVecB=0, gVecX=0;
convert<MatrixType>(symm, gSymm);
convert<MatrixType>(symm, gMatA);
convert<VectorType>(vecB, gVecB);
convert<VectorType>(vecB, gVecX);
Gsl::cholesky(gMatA);
Gsl::cholesky_solve(gMatA, gVecB, gVecX);
VectorType vecX(rows), _vecX, _vecB;
convert(gVecX, _vecX);
symm.llt().solve(vecB, &vecX);
Gsl::prod(gSymm, gVecX, gVecB);
convert(gVecB, _vecB);
// test gsl itself !
VERIFY_IS_APPROX(vecB, _vecB);
VERIFY_IS_APPROX(vecX, _vecX);
Gsl::free(gMatA);
Gsl::free(gSymm);
Gsl::free(gVecB);
Gsl::free(gVecX);
}
#endif
{
LLT<SquareMatrixType> chol(symm);
VERIFY(chol.isPositiveDefinite());
VERIFY_IS_APPROX(symm, chol.matrixL() * chol.matrixL().adjoint());
chol.solve(vecB, &vecX);
VERIFY_IS_APPROX(symm * vecX, vecB);
chol.solve(matB, &matX);
VERIFY_IS_APPROX(symm * matX, matB);
}
int sign = ei_random<int>()%2 ? 1 : -1;
if(sign == -1)
{
symm = -symm; // test a negative matrix
}
{
LDLT<SquareMatrixType> ldlt(symm);
VERIFY(ldlt.isInvertible());
if(sign == 1)
{
VERIFY(ldlt.isPositive());
VERIFY(ldlt.isPositiveDefinite());
}
if(sign == -1)
{
VERIFY(ldlt.isNegative());
VERIFY(ldlt.isNegativeDefinite());
}
// TODO(keir): This doesn't make sense now that LDLT pivots.
//VERIFY_IS_APPROX(symm, ldlt.matrixL() * ldlt.vectorD().asDiagonal() * ldlt.matrixL().adjoint());
ldlt.solve(vecB, &vecX);
VERIFY_IS_APPROX(symm * vecX, vecB);
ldlt.solve(matB, &matX);
VERIFY_IS_APPROX(symm * matX, matB);
}
// test isPositiveDefinite on non definite matrix
if (rows>4)
{
SquareMatrixType symm = a0.block(0,0,rows,cols-4) * a0.block(0,0,rows,cols-4).adjoint();
LLT<SquareMatrixType> chol(symm);
VERIFY(!chol.isPositiveDefinite());
LDLT<SquareMatrixType> cholnosqrt(symm);
VERIFY(!cholnosqrt.isPositiveDefinite());
}
}
template<typename Derived>
void doSomeRankPreservingOperations(Eigen::MatrixBase<Derived>& m)
{
typedef typename Derived::RealScalar RealScalar;
for(int a = 0; a < 3*(m.rows()+m.cols()); a++)
{
RealScalar d = Eigen::ei_random<RealScalar>(-1,1);
int i = Eigen::ei_random<int>(0,m.rows()-1); // i is a random row number
int j;
do {
j = Eigen::ei_random<int>(0,m.rows()-1);
} while (i==j); // j is another one (must be different)
m.row(i) += d * m.row(j);
i = Eigen::ei_random<int>(0,m.cols()-1); // i is a random column number
do {
j = Eigen::ei_random<int>(0,m.cols()-1);
} while (i==j); // j is another one (must be different)
m.col(i) += d * m.col(j);
}
}
template<typename MatrixType> void ldlt_rank()
{
// NOTE there seems to be a problem with too small sizes -- could easily lie in the doSomeRankPreservingOperations function
int rows = ei_random<int>(50,200);
int rank = ei_random<int>(40, rows-1);
// generate a random positive matrix a of given rank
MatrixType m = MatrixType::Random(rows,rows);
QR<MatrixType> qr(m);
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> DiagVectorType;
DiagVectorType d(rows);
d.setZero();
for(int i = 0; i < rank; i++) d(i)=RealScalar(1);
MatrixType a = qr.matrixQ() * d.asDiagonal() * qr.matrixQ().adjoint();
LDLT<MatrixType> ldlt(a);
VERIFY( ei_abs(ldlt.rank() - rank) <= rank / 20 ); // yes, LDLT::rank is a bit inaccurate...
}
void test_cholesky()
{
for(int i = 0; i < g_repeat; i++) {
CALL_SUBTEST( cholesky(Matrix<double,1,1>()) );
CALL_SUBTEST( cholesky(Matrix2d()) );
CALL_SUBTEST( cholesky(Matrix3f()) );
CALL_SUBTEST( cholesky(Matrix4d()) );
CALL_SUBTEST( cholesky(MatrixXcd(7,7)) );
CALL_SUBTEST( cholesky(MatrixXd(17,17)) );
CALL_SUBTEST( cholesky(MatrixXf(200,200)) );
}
for(int i = 0; i < g_repeat/3; i++) {
CALL_SUBTEST( ldlt_rank<MatrixXd>() );
CALL_SUBTEST( ldlt_rank<MatrixXf>() );
CALL_SUBTEST( ldlt_rank<MatrixXcd>() );
}
}