// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2008-2009 Benoit Jacob // // This Source Code Form is subject to the terms of the Mozilla // Public License v. 2.0. If a copy of the MPL was not distributed // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. #include "main.h" #include #include "solverbase.h" using namespace std; template typename MatrixType::RealScalar matrix_l1_norm(const MatrixType& m) { return m.cwiseAbs().colwise().sum().maxCoeff(); } template void lu_non_invertible() { typedef typename MatrixType::RealScalar RealScalar; /* this test covers the following files: LU.h */ Index rows, cols, cols2; if (MatrixType::RowsAtCompileTime == Dynamic) { rows = internal::random(2, EIGEN_TEST_MAX_SIZE); } else { rows = MatrixType::RowsAtCompileTime; } if (MatrixType::ColsAtCompileTime == Dynamic) { cols = internal::random(2, EIGEN_TEST_MAX_SIZE); cols2 = internal::random(2, EIGEN_TEST_MAX_SIZE); } else { cols2 = cols = MatrixType::ColsAtCompileTime; } enum { RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime }; typedef typename internal::kernel_retval_base >::ReturnType KernelMatrixType; typedef typename internal::image_retval_base >::ReturnType ImageMatrixType; typedef Matrix CMatrixType; typedef Matrix RMatrixType; Index rank = internal::random(1, (std::min)(rows, cols) - 1); // The image of the zero matrix should consist of a single (zero) column vector VERIFY((MatrixType::Zero(rows, cols).fullPivLu().image(MatrixType::Zero(rows, cols)).cols() == 1)); // The kernel of the zero matrix is the entire space, and thus is an invertible matrix of dimensions cols. KernelMatrixType kernel = MatrixType::Zero(rows, cols).fullPivLu().kernel(); VERIFY((kernel.fullPivLu().isInvertible())); MatrixType m1(rows, cols), m3(rows, cols2); CMatrixType m2(cols, cols2); createRandomPIMatrixOfRank(rank, rows, cols, m1); FullPivLU lu; // The special value 0.01 below works well in tests. Keep in mind that we're only computing the rank // of singular values are either 0 or 1. // So it's not clear at all that the epsilon should play any role there. lu.setThreshold(RealScalar(0.01)); lu.compute(m1); MatrixType u(rows, cols); u = lu.matrixLU().template triangularView(); RMatrixType l = RMatrixType::Identity(rows, rows); l.block(0, 0, rows, (std::min)(rows, cols)).template triangularView() = lu.matrixLU().block(0, 0, rows, (std::min)(rows, cols)); VERIFY_IS_APPROX(lu.permutationP() * m1 * lu.permutationQ(), l * u); KernelMatrixType m1kernel = lu.kernel(); ImageMatrixType m1image = lu.image(m1); VERIFY_IS_APPROX(m1, lu.reconstructedMatrix()); VERIFY(rank == lu.rank()); VERIFY(cols - lu.rank() == lu.dimensionOfKernel()); VERIFY(!lu.isInjective()); VERIFY(!lu.isInvertible()); VERIFY(!lu.isSurjective()); VERIFY_IS_MUCH_SMALLER_THAN((m1 * m1kernel), m1); VERIFY(m1image.fullPivLu().rank() == rank); VERIFY_IS_APPROX(m1 * m1.adjoint() * m1image, m1image); check_solverbase(m1, lu, rows, cols, cols2); m2 = CMatrixType::Random(cols, cols2); m3 = m1 * m2; m2 = CMatrixType::Random(cols, cols2); // test that the code, which does resize(), may be applied to an xpr m2.block(0, 0, m2.rows(), m2.cols()) = lu.solve(m3); VERIFY_IS_APPROX(m3, m1 * m2); } template void lu_invertible() { /* this test covers the following files: FullPivLU.h */ typedef typename NumTraits::Real RealScalar; Index size = MatrixType::RowsAtCompileTime; if (size == Dynamic) size = internal::random(1, EIGEN_TEST_MAX_SIZE); MatrixType m1(size, size), m2(size, size), m3(size, size); FullPivLU lu; lu.setThreshold(RealScalar(0.01)); do { m1 = MatrixType::Random(size, size); lu.compute(m1); } while (!lu.isInvertible()); VERIFY_IS_APPROX(m1, lu.reconstructedMatrix()); VERIFY(0 == lu.dimensionOfKernel()); VERIFY(lu.kernel().cols() == 1); // the kernel() should consist of a single (zero) column vector VERIFY(size == lu.rank()); VERIFY(lu.isInjective()); VERIFY(lu.isSurjective()); VERIFY(lu.isInvertible()); VERIFY(lu.image(m1).fullPivLu().isInvertible()); check_solverbase(m1, lu, size, size, size); MatrixType m1_inverse = lu.inverse(); m3 = MatrixType::Random(size, size); m2 = lu.solve(m3); VERIFY_IS_APPROX(m2, m1_inverse * m3); RealScalar rcond = (RealScalar(1) / matrix_l1_norm(m1)) / matrix_l1_norm(m1_inverse); const RealScalar rcond_est = lu.rcond(); // Verify that the estimated condition number is within a factor of 10 of the // truth. VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10); // Regression test for Bug 302 MatrixType m4 = MatrixType::Random(size, size); VERIFY_IS_APPROX(lu.solve(m3 * m4), lu.solve(m3) * m4); } template void lu_partial_piv(Index size = MatrixType::ColsAtCompileTime) { /* this test covers the following files: PartialPivLU.h */ typedef typename NumTraits::Real RealScalar; MatrixType m1(size, size), m2(size, size), m3(size, size); m1.setRandom(); PartialPivLU plu(m1); VERIFY_IS_APPROX(m1, plu.reconstructedMatrix()); check_solverbase(m1, plu, size, size, size); MatrixType m1_inverse = plu.inverse(); m3 = MatrixType::Random(size, size); m2 = plu.solve(m3); VERIFY_IS_APPROX(m2, m1_inverse * m3); RealScalar rcond = (RealScalar(1) / matrix_l1_norm(m1)) / matrix_l1_norm(m1_inverse); const RealScalar rcond_est = plu.rcond(); // Verify that the estimate is within a factor of 10 of the truth. VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10); } template void lu_verify_assert() { MatrixType tmp; FullPivLU lu; VERIFY_RAISES_ASSERT(lu.matrixLU()) VERIFY_RAISES_ASSERT(lu.permutationP()) VERIFY_RAISES_ASSERT(lu.permutationQ()) VERIFY_RAISES_ASSERT(lu.kernel()) VERIFY_RAISES_ASSERT(lu.image(tmp)) VERIFY_RAISES_ASSERT(lu.solve(tmp)) VERIFY_RAISES_ASSERT(lu.transpose().solve(tmp)) VERIFY_RAISES_ASSERT(lu.adjoint().solve(tmp)) VERIFY_RAISES_ASSERT(lu.determinant()) VERIFY_RAISES_ASSERT(lu.rank()) VERIFY_RAISES_ASSERT(lu.dimensionOfKernel()) VERIFY_RAISES_ASSERT(lu.isInjective()) VERIFY_RAISES_ASSERT(lu.isSurjective()) VERIFY_RAISES_ASSERT(lu.isInvertible()) VERIFY_RAISES_ASSERT(lu.inverse()) PartialPivLU plu; VERIFY_RAISES_ASSERT(plu.matrixLU()) VERIFY_RAISES_ASSERT(plu.permutationP()) VERIFY_RAISES_ASSERT(plu.solve(tmp)) VERIFY_RAISES_ASSERT(plu.transpose().solve(tmp)) VERIFY_RAISES_ASSERT(plu.adjoint().solve(tmp)) VERIFY_RAISES_ASSERT(plu.determinant()) VERIFY_RAISES_ASSERT(plu.inverse()) } // Rank-deficient matrix returns 0. // https://gitlab.com/libeigen/eigen/-/issues/2889 void test_2889() { Eigen::MatrixXd A = Eigen::MatrixXd({{0.0000000000000000, 0.0000000000000000, 1.0000000000000000, 0.0000000000000000, 0.34149999916553497, 0.0000000000000000, 0.79877008515664061}, {0.0000000000000000, 1.0000000000000000, 0.0000000000000000, 0.29200000315904617, 0.0000000000000000, -0.37149999849498272, -0.16425902650844920}, {0.0000000000000000, 0.0000000000000000, 1.0000000000000000, 0.0000000000000000, 0.34149999916553497, 0.0000000000000000, 0.79877008515664061}, {0.0000000000000000, 1.0000000000000000, 0.0000000000000000, 0.040500000119209290, 0.0000000000000000, -0.30099999904632568, -0.081170580429391403}, {1.0000000000000000, 0.0000000000000000, 0.0000000000000000, 0.0000000000000000, 0.0000000000000000, 0.0000000000000000, -0.0000000000000000}, {0.0000000000000000, 0.70710689672598170, 0.70710666564709435, 0.027000000700354562, 0.025455838867477515, -0.025455847186317101, -0.0068972271572272821}, {1.0000000000000000, 0.0000000000000000, 0.0000000000000000, 0.0000000000000000, 0.0000000000000000, 0.0000000000000000, -0.0000000000000000}}); Eigen::FullPivLU lu_factorization(A); double rcond = lu_factorization.rcond(); VERIFY_IS_EQUAL(rcond, 0.0); } EIGEN_DECLARE_TEST(lu) { for (int i = 0; i < g_repeat; i++) { CALL_SUBTEST_1(lu_non_invertible()); CALL_SUBTEST_1(lu_invertible()); CALL_SUBTEST_1(lu_verify_assert()); CALL_SUBTEST_1(lu_partial_piv()); CALL_SUBTEST_2((lu_non_invertible >())); CALL_SUBTEST_2((lu_verify_assert >())); CALL_SUBTEST_2(lu_partial_piv()); CALL_SUBTEST_2(lu_partial_piv()); CALL_SUBTEST_2((lu_partial_piv >())); CALL_SUBTEST_3(lu_non_invertible()); CALL_SUBTEST_3(lu_invertible()); CALL_SUBTEST_3(lu_verify_assert()); CALL_SUBTEST_4(lu_non_invertible()); CALL_SUBTEST_4(lu_invertible()); CALL_SUBTEST_4(lu_partial_piv(internal::random(1, EIGEN_TEST_MAX_SIZE))); CALL_SUBTEST_4(lu_verify_assert()); CALL_SUBTEST_5(lu_non_invertible()); CALL_SUBTEST_5(lu_invertible()); CALL_SUBTEST_5(lu_verify_assert()); CALL_SUBTEST_6(lu_non_invertible()); CALL_SUBTEST_6(lu_invertible()); CALL_SUBTEST_6(lu_partial_piv(internal::random(1, EIGEN_TEST_MAX_SIZE))); CALL_SUBTEST_6(lu_verify_assert()); CALL_SUBTEST_7((lu_non_invertible >())); // Test problem size constructors CALL_SUBTEST_8(PartialPivLU(10)); CALL_SUBTEST_8(FullPivLU(10, 20);); CALL_SUBTEST_9(test_2889()); } }