From 80b513378948f78bc7729c431eb68ca5513a1d62 Mon Sep 17 00:00:00 2001 From: Gael Guennebaud Date: Thu, 6 Oct 2016 09:55:50 +0200 Subject: [PATCH] Fix compilation of qr.inverse() for column and full pivoting variants. --- Eigen/src/QR/ColPivHouseholderQR.h | 6 +++--- Eigen/src/QR/FullPivHouseholderQR.h | 6 +++--- test/qr_colpivoting.cpp | 12 ++++++++++++ test/qr_fullpivoting.cpp | 12 ++++++++++++ 4 files changed, 30 insertions(+), 6 deletions(-) diff --git a/Eigen/src/QR/ColPivHouseholderQR.h b/Eigen/src/QR/ColPivHouseholderQR.h index 9650781d6..618f80480 100644 --- a/Eigen/src/QR/ColPivHouseholderQR.h +++ b/Eigen/src/QR/ColPivHouseholderQR.h @@ -613,12 +613,12 @@ void ColPivHouseholderQR<_MatrixType>::_solve_impl(const RhsType &rhs, DstType & namespace internal { -template -struct Assignment >, internal::assign_op, Dense2Dense> +template +struct Assignment >, internal::assign_op::Scalar>, Dense2Dense> { typedef ColPivHouseholderQR QrType; typedef Inverse SrcXprType; - static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op &) + static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op &) { dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols())); } diff --git a/Eigen/src/QR/FullPivHouseholderQR.h b/Eigen/src/QR/FullPivHouseholderQR.h index e0e15100d..e489bddc2 100644 --- a/Eigen/src/QR/FullPivHouseholderQR.h +++ b/Eigen/src/QR/FullPivHouseholderQR.h @@ -575,12 +575,12 @@ void FullPivHouseholderQR<_MatrixType>::_solve_impl(const RhsType &rhs, DstType namespace internal { -template -struct Assignment >, internal::assign_op, Dense2Dense> +template +struct Assignment >, internal::assign_op::Scalar>, Dense2Dense> { typedef FullPivHouseholderQR QrType; typedef Inverse SrcXprType; - static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op &) + static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op &) { dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols())); } diff --git a/test/qr_colpivoting.cpp b/test/qr_colpivoting.cpp index 057bb014c..26ed27f5c 100644 --- a/test/qr_colpivoting.cpp +++ b/test/qr_colpivoting.cpp @@ -141,6 +141,18 @@ template void qr() m2 = MatrixType::Random(cols,cols2); m2 = qr.solve(m3); VERIFY_IS_APPROX(m3, m1*m2); + + { + Index size = rows; + do { + m1 = MatrixType::Random(size,size); + qr.compute(m1); + } while(!qr.isInvertible()); + MatrixType m1_inv = qr.inverse(); + m3 = m1 * MatrixType::Random(size,cols2); + m2 = qr.solve(m3); + VERIFY_IS_APPROX(m2, m1_inv*m3); + } } template void qr_fixedsize() diff --git a/test/qr_fullpivoting.cpp b/test/qr_fullpivoting.cpp index 05a705887..70e89c198 100644 --- a/test/qr_fullpivoting.cpp +++ b/test/qr_fullpivoting.cpp @@ -54,6 +54,18 @@ template void qr() m2 = MatrixType::Random(cols,cols2); m2 = qr.solve(m3); VERIFY_IS_APPROX(m3, m1*m2); + + { + Index size = rows; + do { + m1 = MatrixType::Random(size,size); + qr.compute(m1); + } while(!qr.isInvertible()); + MatrixType m1_inv = qr.inverse(); + m3 = m1 * MatrixType::Random(size,cols2); + m2 = qr.solve(m3); + VERIFY_IS_APPROX(m2, m1_inv*m3); + } } template void qr_invertible()