Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 1 | // This file is part of Eigen, a lightweight C++ template library |
| 2 | // for linear algebra. |
| 3 | // |
| 4 | // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr> |
| 5 | // |
| 6 | // This Source Code Form is subject to the terms of the Mozilla |
| 7 | // Public License v. 2.0. If a copy of the MPL was not distributed |
| 8 | // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. |
| 9 | |
| 10 | #include "main.h" |
| 11 | #include <Eigen/QR> |
| 12 | |
| 13 | template<typename MatrixType> void qr(const MatrixType& m) |
| 14 | { |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 15 | Index rows = m.rows(); |
| 16 | Index cols = m.cols(); |
| 17 | |
| 18 | typedef typename MatrixType::Scalar Scalar; |
| 19 | typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType; |
| 20 | |
| 21 | MatrixType a = MatrixType::Random(rows,cols); |
| 22 | HouseholderQR<MatrixType> qrOfA(a); |
| 23 | |
| 24 | MatrixQType q = qrOfA.householderQ(); |
| 25 | VERIFY_IS_UNITARY(q); |
| 26 | |
| 27 | MatrixType r = qrOfA.matrixQR().template triangularView<Upper>(); |
| 28 | VERIFY_IS_APPROX(a, qrOfA.householderQ() * r); |
| 29 | } |
| 30 | |
| 31 | template<typename MatrixType, int Cols2> void qr_fixedsize() |
| 32 | { |
| 33 | enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime }; |
| 34 | typedef typename MatrixType::Scalar Scalar; |
| 35 | Matrix<Scalar,Rows,Cols> m1 = Matrix<Scalar,Rows,Cols>::Random(); |
| 36 | HouseholderQR<Matrix<Scalar,Rows,Cols> > qr(m1); |
| 37 | |
| 38 | Matrix<Scalar,Rows,Cols> r = qr.matrixQR(); |
| 39 | // FIXME need better way to construct trapezoid |
| 40 | for(int i = 0; i < Rows; i++) for(int j = 0; j < Cols; j++) if(i>j) r(i,j) = Scalar(0); |
| 41 | |
| 42 | VERIFY_IS_APPROX(m1, qr.householderQ() * r); |
| 43 | |
| 44 | Matrix<Scalar,Cols,Cols2> m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2); |
| 45 | Matrix<Scalar,Rows,Cols2> m3 = m1*m2; |
| 46 | m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2); |
| 47 | m2 = qr.solve(m3); |
| 48 | VERIFY_IS_APPROX(m3, m1*m2); |
| 49 | } |
| 50 | |
| 51 | template<typename MatrixType> void qr_invertible() |
| 52 | { |
| 53 | using std::log; |
| 54 | using std::abs; |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 55 | using std::pow; |
| 56 | using std::max; |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 57 | typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar; |
| 58 | typedef typename MatrixType::Scalar Scalar; |
| 59 | |
| 60 | int size = internal::random<int>(10,50); |
| 61 | |
| 62 | MatrixType m1(size, size), m2(size, size), m3(size, size); |
| 63 | m1 = MatrixType::Random(size,size); |
| 64 | |
| 65 | if (internal::is_same<RealScalar,float>::value) |
| 66 | { |
| 67 | // let's build a matrix more stable to inverse |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 68 | MatrixType a = MatrixType::Random(size,size*4); |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 69 | m1 += a * a.adjoint(); |
| 70 | } |
| 71 | |
| 72 | HouseholderQR<MatrixType> qr(m1); |
| 73 | m3 = MatrixType::Random(size,size); |
| 74 | m2 = qr.solve(m3); |
| 75 | VERIFY_IS_APPROX(m3, m1*m2); |
| 76 | |
| 77 | // now construct a matrix with prescribed determinant |
| 78 | m1.setZero(); |
| 79 | for(int i = 0; i < size; i++) m1(i,i) = internal::random<Scalar>(); |
| 80 | RealScalar absdet = abs(m1.diagonal().prod()); |
| 81 | m3 = qr.householderQ(); // get a unitary |
| 82 | m1 = m3 * m1 * m3; |
| 83 | qr.compute(m1); |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 84 | VERIFY_IS_APPROX(log(absdet), qr.logAbsDeterminant()); |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 85 | // This test is tricky if the determinant becomes too small. |
| 86 | // Since we generate random numbers with magnitude rrange [0,1], the average determinant is 0.5^size |
| 87 | VERIFY_IS_MUCH_SMALLER_THAN( abs(absdet-qr.absDeterminant()), numext::maxi(RealScalar(pow(0.5,size)),numext::maxi<RealScalar>(abs(absdet),abs(qr.absDeterminant()))) ); |
| 88 | |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 89 | } |
| 90 | |
| 91 | template<typename MatrixType> void qr_verify_assert() |
| 92 | { |
| 93 | MatrixType tmp; |
| 94 | |
| 95 | HouseholderQR<MatrixType> qr; |
| 96 | VERIFY_RAISES_ASSERT(qr.matrixQR()) |
| 97 | VERIFY_RAISES_ASSERT(qr.solve(tmp)) |
| 98 | VERIFY_RAISES_ASSERT(qr.householderQ()) |
| 99 | VERIFY_RAISES_ASSERT(qr.absDeterminant()) |
| 100 | VERIFY_RAISES_ASSERT(qr.logAbsDeterminant()) |
| 101 | } |
| 102 | |
| 103 | void test_qr() |
| 104 | { |
| 105 | for(int i = 0; i < g_repeat; i++) { |
| 106 | CALL_SUBTEST_1( qr(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
| 107 | CALL_SUBTEST_2( qr(MatrixXcd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2),internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))) ); |
| 108 | CALL_SUBTEST_3(( qr_fixedsize<Matrix<float,3,4>, 2 >() )); |
| 109 | CALL_SUBTEST_4(( qr_fixedsize<Matrix<double,6,2>, 4 >() )); |
| 110 | CALL_SUBTEST_5(( qr_fixedsize<Matrix<double,2,5>, 7 >() )); |
| 111 | CALL_SUBTEST_11( qr(Matrix<float,1,1>()) ); |
| 112 | } |
| 113 | |
| 114 | for(int i = 0; i < g_repeat; i++) { |
| 115 | CALL_SUBTEST_1( qr_invertible<MatrixXf>() ); |
| 116 | CALL_SUBTEST_6( qr_invertible<MatrixXd>() ); |
| 117 | CALL_SUBTEST_7( qr_invertible<MatrixXcf>() ); |
| 118 | CALL_SUBTEST_8( qr_invertible<MatrixXcd>() ); |
| 119 | } |
| 120 | |
| 121 | CALL_SUBTEST_9(qr_verify_assert<Matrix3f>()); |
| 122 | CALL_SUBTEST_10(qr_verify_assert<Matrix3d>()); |
| 123 | CALL_SUBTEST_1(qr_verify_assert<MatrixXf>()); |
| 124 | CALL_SUBTEST_6(qr_verify_assert<MatrixXd>()); |
| 125 | CALL_SUBTEST_7(qr_verify_assert<MatrixXcf>()); |
| 126 | CALL_SUBTEST_8(qr_verify_assert<MatrixXcd>()); |
| 127 | |
| 128 | // Test problem size constructors |
| 129 | CALL_SUBTEST_12(HouseholderQR<MatrixXf>(10, 20)); |
| 130 | } |