aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/qr_colpivoting.cpp
diff options
context:
space:
mode:
authorGravatar Rasmus Munk Larsen <rmlarsen@google.com>2016-01-28 15:07:26 -0800
committerGravatar Rasmus Munk Larsen <rmlarsen@google.com>2016-01-28 15:07:26 -0800
commitacce4dd0500fbb9524fe35aacafb7fbc5f7f76f9 (patch)
tree56a017d84a84c45564278b6c284ec4387963618d /test/qr_colpivoting.cpp
parentb908e071a80fce910efc82c1c50dd6be1e226dcd (diff)
Change Eigen's ColPivHouseholderQR to use the numerically stable norm downdate formula from http://www.netlib.org/lapack/lawnspdf/lawn176.pdf, which has been used in LAPACK's xGEQPF and xGEQP3 since 2006. With the old formula, the code chooses the wrong pivots and fails to correctly determine rank on graded matrices.
This change also adds additional checks for non-increasing diagonal in R11 to existing unit tests, and adds a new unit test with the Kahan matrix, which consistently fails for the original code. Benchmark timings on Intel(R) Xeon(R) CPU E5-1650 v3 @ 3.50GHz. Code compiled with AVX & FMA. I just ran on square matrices of 3 difference sizes. Benchmark Time(ns) CPU(ns) Iterations ------------------------------------------------------- Before: BM_EigencolPivQR/64 53677 53627 12890 BM_EigencolPivQR/512 15265408 15250784 46 BM_EigencolPivQR/4k 15403556228 15388788368 2 After (non-vectorized version): Benchmark Time(ns) CPU(ns) Iterations Degradation -------------------------------------------------------------------- BM_EigencolPivQR/64 63736 63669 10844 18.5% BM_EigencolPivQR/512 16052546 16037381 43 5.1% BM_EigencolPivQR/4k 15149263620 15132025316 2 -2.0% Performance-wise there seems to be a ~18.5% degradation for small (64x64) matrices, probably due to the cost of more O(min(m,n)^2) sqrt operations that are not needed for the unstable formula.
Diffstat (limited to 'test/qr_colpivoting.cpp')
-rw-r--r--test/qr_colpivoting.cpp86
1 files changed, 86 insertions, 0 deletions
diff --git a/test/qr_colpivoting.cpp b/test/qr_colpivoting.cpp
index e7abd3725..7b97292db 100644
--- a/test/qr_colpivoting.cpp
+++ b/test/qr_colpivoting.cpp
@@ -19,6 +19,7 @@ template<typename MatrixType> void qr()
Index rank = internal::random<Index>(1, (std::min)(rows, cols)-1);
typedef typename MatrixType::Scalar Scalar;
+ typedef typename MatrixType::RealScalar RealScalar;
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
MatrixType m1;
createRandomPIMatrixOfRank(rank,rows,cols,m1);
@@ -36,6 +37,24 @@ template<typename MatrixType> void qr()
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<Scalar>::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(x, y)) {
+ 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);
@@ -47,6 +66,7 @@ template<typename MatrixType, int Cols2> void qr_fixedsize()
{
enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime };
typedef typename MatrixType::Scalar Scalar;
+ typedef typename MatrixType::RealScalar RealScalar;
int rank = internal::random<int>(1, (std::min)(int(Rows), int(Cols))-1);
Matrix<Scalar,Rows,Cols> m1;
createRandomPIMatrixOfRank(rank,Rows,Cols,m1);
@@ -66,6 +86,67 @@ template<typename MatrixType, int Cols2> void qr_fixedsize()
m2 = Matrix<Scalar,Cols,Cols2>::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<Scalar>::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(x, y)) {
+ 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<typename MatrixType> 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<RealScalar>::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<MatrixType> qr(m1);
+ MatrixType r = qr.matrixQR().template triangularView<Upper>();
+
+ RealScalar threshold =
+ std::sqrt(RealScalar(rows)) * (std::abs)(r(0, 0)) * NumTraits<Scalar>::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(x, y)) {
+ 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<typename MatrixType> void qr_invertible()
@@ -132,6 +213,11 @@ void test_qr_colpivoting()
}
for(int i = 0; i < g_repeat; i++) {
+ CALL_SUBTEST_1( qr_kahan_matrix<MatrixXf>() );
+ CALL_SUBTEST_2( qr_kahan_matrix<MatrixXd>() );
+ }
+
+ for(int i = 0; i < g_repeat; i++) {
CALL_SUBTEST_1( qr_invertible<MatrixXf>() );
CALL_SUBTEST_2( qr_invertible<MatrixXd>() );
CALL_SUBTEST_6( qr_invertible<MatrixXcf>() );