// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2008 Gael Guennebaud // Copyright (C) 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 template void cod() { typedef typename MatrixType::Index Index; Index rows = internal::random(2, EIGEN_TEST_MAX_SIZE); Index cols = internal::random(2, EIGEN_TEST_MAX_SIZE); Index cols2 = internal::random(2, EIGEN_TEST_MAX_SIZE); Index rank = internal::random(1, (std::min)(rows, cols) - 1); typedef typename MatrixType::Scalar Scalar; typedef Matrix MatrixQType; MatrixType matrix; createRandomPIMatrixOfRank(rank, rows, cols, matrix); CompleteOrthogonalDecomposition cod(matrix); VERIFY(rank == cod.rank()); VERIFY(cols - cod.rank() == cod.dimensionOfKernel()); VERIFY(!cod.isInjective()); VERIFY(!cod.isInvertible()); VERIFY(!cod.isSurjective()); MatrixQType q = cod.householderQ(); VERIFY_IS_UNITARY(q); MatrixType z = cod.matrixZ(); VERIFY_IS_UNITARY(z); MatrixType t; t.setZero(rows, cols); t.topLeftCorner(rank, rank) = cod.matrixT().topLeftCorner(rank, rank).template triangularView(); MatrixType c = q * t * z * cod.colsPermutation().inverse(); VERIFY_IS_APPROX(matrix, c); MatrixType exact_solution = MatrixType::Random(cols, cols2); MatrixType rhs = matrix * exact_solution; MatrixType cod_solution = cod.solve(rhs); VERIFY_IS_APPROX(rhs, matrix * cod_solution); // Verify that we get the same minimum-norm solution as the SVD. JacobiSVD svd(matrix, ComputeThinU | ComputeThinV); MatrixType svd_solution = svd.solve(rhs); VERIFY_IS_APPROX(cod_solution, svd_solution); } template void cod_fixedsize() { enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime }; typedef typename MatrixType::Scalar Scalar; int rank = internal::random(1, (std::min)(int(Rows), int(Cols)) - 1); Matrix matrix; createRandomPIMatrixOfRank(rank, Rows, Cols, matrix); CompleteOrthogonalDecomposition > cod(matrix); VERIFY(rank == cod.rank()); VERIFY(Cols - cod.rank() == cod.dimensionOfKernel()); VERIFY(cod.isInjective() == (rank == Rows)); VERIFY(cod.isSurjective() == (rank == Cols)); VERIFY(cod.isInvertible() == (cod.isInjective() && cod.isSurjective())); Matrix exact_solution; exact_solution.setRandom(Cols, Cols2); Matrix rhs = matrix * exact_solution; Matrix cod_solution = cod.solve(rhs); VERIFY_IS_APPROX(rhs, matrix * cod_solution); // Verify that we get the same minimum-norm solution as the SVD. JacobiSVD svd(matrix, ComputeFullU | ComputeFullV); Matrix svd_solution = svd.solve(rhs); VERIFY_IS_APPROX(cod_solution, svd_solution); } template void qr() { typedef typename MatrixType::Index Index; Index rows = internal::random(2,EIGEN_TEST_MAX_SIZE), cols = internal::random(2,EIGEN_TEST_MAX_SIZE), cols2 = internal::random(2,EIGEN_TEST_MAX_SIZE); Index rank = internal::random(1, (std::min)(rows, cols)-1); typedef typename MatrixType::Scalar Scalar; typedef typename MatrixType::RealScalar RealScalar; typedef Matrix MatrixQType; MatrixType m1; createRandomPIMatrixOfRank(rank,rows,cols,m1); ColPivHouseholderQR qr(m1); VERIFY_IS_EQUAL(rank, qr.rank()); VERIFY_IS_EQUAL(cols - qr.rank(), qr.dimensionOfKernel()); VERIFY(!qr.isInjective()); VERIFY(!qr.isInvertible()); VERIFY(!qr.isSurjective()); MatrixQType q = qr.householderQ(); VERIFY_IS_UNITARY(q); MatrixType r = qr.matrixQR().template triangularView(); MatrixType c = q * r * qr.colsPermutation().inverse(); VERIFY_IS_APPROX(m1, c); // Verify that the absolute value of the diagonal elements in R are // non-increasing until they reach the singularity threshold. RealScalar threshold = std::sqrt(RealScalar(rows)) * (std::abs)(r(0, 0)) * NumTraits::epsilon(); for (Index i = 0; i < (std::min)(rows, cols) - 1; ++i) { RealScalar x = (std::abs)(r(i, i)); RealScalar y = (std::abs)(r(i + 1, i + 1)); if (x < threshold && y < threshold) continue; if (!test_isApproxOrLessThan(y, x)) { for (Index j = 0; j < (std::min)(rows, cols); ++j) { std::cout << "i = " << j << ", |r_ii| = " << (std::abs)(r(j, j)) << std::endl; } std::cout << "Failure at i=" << i << ", rank=" << rank << ", threshold=" << threshold << std::endl; } VERIFY_IS_APPROX_OR_LESS_THAN(y, x); } MatrixType m2 = MatrixType::Random(cols,cols2); MatrixType m3 = m1*m2; m2 = MatrixType::Random(cols,cols2); m2 = qr.solve(m3); VERIFY_IS_APPROX(m3, m1*m2); } template void qr_fixedsize() { enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime }; typedef typename MatrixType::Scalar Scalar; typedef typename MatrixType::RealScalar RealScalar; int rank = internal::random(1, (std::min)(int(Rows), int(Cols))-1); Matrix m1; createRandomPIMatrixOfRank(rank,Rows,Cols,m1); ColPivHouseholderQR > qr(m1); VERIFY_IS_EQUAL(rank, qr.rank()); VERIFY_IS_EQUAL(Cols - qr.rank(), qr.dimensionOfKernel()); VERIFY_IS_EQUAL(qr.isInjective(), (rank == Rows)); VERIFY_IS_EQUAL(qr.isSurjective(), (rank == Cols)); VERIFY_IS_EQUAL(qr.isInvertible(), (qr.isInjective() && qr.isSurjective())); Matrix r = qr.matrixQR().template triangularView(); Matrix c = qr.householderQ() * r * qr.colsPermutation().inverse(); VERIFY_IS_APPROX(m1, c); Matrix m2 = Matrix::Random(Cols,Cols2); Matrix m3 = m1*m2; m2 = Matrix::Random(Cols,Cols2); m2 = qr.solve(m3); VERIFY_IS_APPROX(m3, m1*m2); // Verify that the absolute value of the diagonal elements in R are // non-increasing until they reache the singularity threshold. RealScalar threshold = std::sqrt(RealScalar(Rows)) * (std::abs)(r(0, 0)) * NumTraits::epsilon(); for (Index i = 0; i < (std::min)(int(Rows), int(Cols)) - 1; ++i) { RealScalar x = (std::abs)(r(i, i)); RealScalar y = (std::abs)(r(i + 1, i + 1)); if (x < threshold && y < threshold) continue; if (!test_isApproxOrLessThan(y, x)) { for (Index j = 0; j < (std::min)(int(Rows), int(Cols)); ++j) { std::cout << "i = " << j << ", |r_ii| = " << (std::abs)(r(j, j)) << std::endl; } std::cout << "Failure at i=" << i << ", rank=" << rank << ", threshold=" << threshold << std::endl; } VERIFY_IS_APPROX_OR_LESS_THAN(y, x); } } // This test is meant to verify that pivots are chosen such that // even for a graded matrix, the diagonal of R falls of roughly // monotonically until it reaches the threshold for singularity. // We use the so-called Kahan matrix, which is a famous counter-example // for rank-revealing QR. See // http://www.netlib.org/lapack/lawnspdf/lawn176.pdf // page 3 for more detail. template void qr_kahan_matrix() { typedef typename MatrixType::Index Index; typedef typename MatrixType::Scalar Scalar; typedef typename MatrixType::RealScalar RealScalar; Index rows = 300, cols = rows; MatrixType m1; m1.setZero(rows,cols); RealScalar s = std::pow(NumTraits::epsilon(), 1.0 / rows); RealScalar c = std::sqrt(1 - s*s); for (Index i = 0; i < rows; ++i) { m1(i, i) = pow(s, i); m1.row(i).tail(rows - i - 1) = -pow(s, i) * c * MatrixType::Ones(1, rows - i - 1); } m1 = (m1 + m1.transpose()).eval(); ColPivHouseholderQR qr(m1); MatrixType r = qr.matrixQR().template triangularView(); RealScalar threshold = std::sqrt(RealScalar(rows)) * (std::abs)(r(0, 0)) * NumTraits::epsilon(); for (Index i = 0; i < (std::min)(rows, cols) - 1; ++i) { RealScalar x = (std::abs)(r(i, i)); RealScalar y = (std::abs)(r(i + 1, i + 1)); if (x < threshold && y < threshold) continue; if (!test_isApproxOrLessThan(y, x)) { for (Index j = 0; j < (std::min)(rows, cols); ++j) { std::cout << "i = " << j << ", |r_ii| = " << (std::abs)(r(j, j)) << std::endl; } std::cout << "Failure at i=" << i << ", rank=" << qr.rank() << ", threshold=" << threshold << std::endl; } VERIFY_IS_APPROX_OR_LESS_THAN(y, x); } } template void qr_invertible() { using std::log; using std::abs; typedef typename NumTraits::Real RealScalar; typedef typename MatrixType::Scalar Scalar; int size = internal::random(10,50); MatrixType m1(size, size), m2(size, size), m3(size, size); m1 = MatrixType::Random(size,size); if (internal::is_same::value) { // let's build a matrix more stable to inverse MatrixType a = MatrixType::Random(size,size*2); m1 += a * a.adjoint(); } ColPivHouseholderQR qr(m1); m3 = MatrixType::Random(size,size); m2 = qr.solve(m3); //VERIFY_IS_APPROX(m3, m1*m2); // now construct a matrix with prescribed determinant m1.setZero(); for(int i = 0; i < size; i++) m1(i,i) = internal::random(); RealScalar absdet = abs(m1.diagonal().prod()); m3 = qr.householderQ(); // get a unitary m1 = m3 * m1 * m3; qr.compute(m1); VERIFY_IS_APPROX(absdet, qr.absDeterminant()); VERIFY_IS_APPROX(log(absdet), qr.logAbsDeterminant()); } template void qr_verify_assert() { MatrixType tmp; ColPivHouseholderQR qr; VERIFY_RAISES_ASSERT(qr.matrixQR()) VERIFY_RAISES_ASSERT(qr.solve(tmp)) VERIFY_RAISES_ASSERT(qr.householderQ()) VERIFY_RAISES_ASSERT(qr.dimensionOfKernel()) VERIFY_RAISES_ASSERT(qr.isInjective()) VERIFY_RAISES_ASSERT(qr.isSurjective()) VERIFY_RAISES_ASSERT(qr.isInvertible()) VERIFY_RAISES_ASSERT(qr.inverse()) VERIFY_RAISES_ASSERT(qr.absDeterminant()) VERIFY_RAISES_ASSERT(qr.logAbsDeterminant()) } void test_qr_colpivoting() { for(int i = 0; i < g_repeat; i++) { CALL_SUBTEST_1( qr() ); CALL_SUBTEST_2( qr() ); CALL_SUBTEST_3( qr() ); CALL_SUBTEST_4(( qr_fixedsize, 4 >() )); CALL_SUBTEST_5(( qr_fixedsize, 3 >() )); CALL_SUBTEST_5(( qr_fixedsize, 1 >() )); } for(int i = 0; i < g_repeat; i++) { CALL_SUBTEST_1( cod() ); CALL_SUBTEST_2( cod() ); CALL_SUBTEST_3( cod() ); CALL_SUBTEST_4(( cod_fixedsize, 4 >() )); CALL_SUBTEST_5(( cod_fixedsize, 3 >() )); CALL_SUBTEST_5(( cod_fixedsize, 1 >() )); } for(int i = 0; i < g_repeat; i++) { CALL_SUBTEST_1( qr_invertible() ); CALL_SUBTEST_2( qr_invertible() ); CALL_SUBTEST_6( qr_invertible() ); CALL_SUBTEST_3( qr_invertible() ); } CALL_SUBTEST_7(qr_verify_assert()); CALL_SUBTEST_8(qr_verify_assert()); CALL_SUBTEST_1(qr_verify_assert()); CALL_SUBTEST_2(qr_verify_assert()); CALL_SUBTEST_6(qr_verify_assert()); CALL_SUBTEST_3(qr_verify_assert()); // Test problem size constructors CALL_SUBTEST_9(ColPivHouseholderQR(10, 20)); CALL_SUBTEST_1( qr_kahan_matrix() ); CALL_SUBTEST_2( qr_kahan_matrix() ); }