// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2008 Gael Guennebaud // // 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 . #define EIGEN_NO_ASSERTION_CHECKING #include "main.h" #include #include #ifdef HAS_GSL #include "gsl_helper.h" #endif template 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::Real RealScalar; typedef Matrix SquareMatrixType; typedef Matrix 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().setZero(); #ifdef HAS_GSL if (ei_is_same_type::ret) { typedef GslTraits Gsl; typename Gsl::Matrix gMatA=0, gSymm=0; typename Gsl::Vector gVecB=0, gVecX=0; convert(symm, gSymm); convert(symm, gMatA); convert(vecB, gVecB); convert(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 chol(symm); 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()%2 ? 1 : -1; if(sign == -1) { symm = -symm; // test a negative matrix } { LDLT ldlt(symm); // 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); } } template void cholesky_verify_assert() { MatrixType tmp; LLT llt; VERIFY_RAISES_ASSERT(llt.matrixL()) VERIFY_RAISES_ASSERT(llt.solve(tmp,&tmp)) VERIFY_RAISES_ASSERT(llt.solveInPlace(&tmp)) LDLT ldlt; VERIFY_RAISES_ASSERT(ldlt.matrixL()) VERIFY_RAISES_ASSERT(ldlt.permutationP()) VERIFY_RAISES_ASSERT(ldlt.vectorD()) VERIFY_RAISES_ASSERT(ldlt.isPositive()) VERIFY_RAISES_ASSERT(ldlt.isNegative()) VERIFY_RAISES_ASSERT(ldlt.solve(tmp,&tmp)) VERIFY_RAISES_ASSERT(ldlt.solveInPlace(&tmp)) } void test_cholesky() { for(int i = 0; i < g_repeat; i++) { CALL_SUBTEST( cholesky(Matrix()) ); CALL_SUBTEST( cholesky(MatrixXd(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)) ); } CALL_SUBTEST( cholesky_verify_assert() ); CALL_SUBTEST( cholesky_verify_assert() ); CALL_SUBTEST( cholesky_verify_assert() ); CALL_SUBTEST( cholesky_verify_assert() ); }