| 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 | {
 | 
|---|
| 15 |   typedef typename MatrixType::Index Index;
 | 
|---|
| 16 | 
 | 
|---|
| 17 |   Index rows = m.rows();
 | 
|---|
| 18 |   Index cols = m.cols();
 | 
|---|
| 19 | 
 | 
|---|
| 20 |   typedef typename MatrixType::Scalar Scalar;
 | 
|---|
| 21 |   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
 | 
|---|
| 22 | 
 | 
|---|
| 23 |   MatrixType a = MatrixType::Random(rows,cols);
 | 
|---|
| 24 |   HouseholderQR<MatrixType> qrOfA(a);
 | 
|---|
| 25 | 
 | 
|---|
| 26 |   MatrixQType q = qrOfA.householderQ();
 | 
|---|
| 27 |   VERIFY_IS_UNITARY(q);
 | 
|---|
| 28 | 
 | 
|---|
| 29 |   MatrixType r = qrOfA.matrixQR().template triangularView<Upper>();
 | 
|---|
| 30 |   VERIFY_IS_APPROX(a, qrOfA.householderQ() * r);
 | 
|---|
| 31 | }
 | 
|---|
| 32 | 
 | 
|---|
| 33 | template<typename MatrixType, int Cols2> void qr_fixedsize()
 | 
|---|
| 34 | {
 | 
|---|
| 35 |   enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime };
 | 
|---|
| 36 |   typedef typename MatrixType::Scalar Scalar;
 | 
|---|
| 37 |   Matrix<Scalar,Rows,Cols> m1 = Matrix<Scalar,Rows,Cols>::Random();
 | 
|---|
| 38 |   HouseholderQR<Matrix<Scalar,Rows,Cols> > qr(m1);
 | 
|---|
| 39 | 
 | 
|---|
| 40 |   Matrix<Scalar,Rows,Cols> r = qr.matrixQR();
 | 
|---|
| 41 |   // FIXME need better way to construct trapezoid
 | 
|---|
| 42 |   for(int i = 0; i < Rows; i++) for(int j = 0; j < Cols; j++) if(i>j) r(i,j) = Scalar(0);
 | 
|---|
| 43 | 
 | 
|---|
| 44 |   VERIFY_IS_APPROX(m1, qr.householderQ() * r);
 | 
|---|
| 45 | 
 | 
|---|
| 46 |   Matrix<Scalar,Cols,Cols2> m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2);
 | 
|---|
| 47 |   Matrix<Scalar,Rows,Cols2> m3 = m1*m2;
 | 
|---|
| 48 |   m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2);
 | 
|---|
| 49 |   m2 = qr.solve(m3);
 | 
|---|
| 50 |   VERIFY_IS_APPROX(m3, m1*m2);
 | 
|---|
| 51 | }
 | 
|---|
| 52 | 
 | 
|---|
| 53 | template<typename MatrixType> void qr_invertible()
 | 
|---|
| 54 | {
 | 
|---|
| 55 |   using std::log;
 | 
|---|
| 56 |   using std::abs;
 | 
|---|
| 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
 | 
|---|
| 68 |     MatrixType a = MatrixType::Random(size,size*2);
 | 
|---|
| 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);
 | 
|---|
| 84 |   VERIFY_IS_APPROX(absdet, qr.absDeterminant());
 | 
|---|
| 85 |   VERIFY_IS_APPROX(log(absdet), qr.logAbsDeterminant());
 | 
|---|
| 86 | }
 | 
|---|
| 87 | 
 | 
|---|
| 88 | template<typename MatrixType> void qr_verify_assert()
 | 
|---|
| 89 | {
 | 
|---|
| 90 |   MatrixType tmp;
 | 
|---|
| 91 | 
 | 
|---|
| 92 |   HouseholderQR<MatrixType> qr;
 | 
|---|
| 93 |   VERIFY_RAISES_ASSERT(qr.matrixQR())
 | 
|---|
| 94 |   VERIFY_RAISES_ASSERT(qr.solve(tmp))
 | 
|---|
| 95 |   VERIFY_RAISES_ASSERT(qr.householderQ())
 | 
|---|
| 96 |   VERIFY_RAISES_ASSERT(qr.absDeterminant())
 | 
|---|
| 97 |   VERIFY_RAISES_ASSERT(qr.logAbsDeterminant())
 | 
|---|
| 98 | }
 | 
|---|
| 99 | 
 | 
|---|
| 100 | void test_qr()
 | 
|---|
| 101 | {
 | 
|---|
| 102 |   for(int i = 0; i < g_repeat; i++) {
 | 
|---|
| 103 |    CALL_SUBTEST_1( qr(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
 | 
|---|
| 104 |    CALL_SUBTEST_2( qr(MatrixXcd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2),internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))) );
 | 
|---|
| 105 |    CALL_SUBTEST_3(( qr_fixedsize<Matrix<float,3,4>, 2 >() ));
 | 
|---|
| 106 |    CALL_SUBTEST_4(( qr_fixedsize<Matrix<double,6,2>, 4 >() ));
 | 
|---|
| 107 |    CALL_SUBTEST_5(( qr_fixedsize<Matrix<double,2,5>, 7 >() ));
 | 
|---|
| 108 |    CALL_SUBTEST_11( qr(Matrix<float,1,1>()) );
 | 
|---|
| 109 |   }
 | 
|---|
| 110 | 
 | 
|---|
| 111 |   for(int i = 0; i < g_repeat; i++) {
 | 
|---|
| 112 |     CALL_SUBTEST_1( qr_invertible<MatrixXf>() );
 | 
|---|
| 113 |     CALL_SUBTEST_6( qr_invertible<MatrixXd>() );
 | 
|---|
| 114 |     CALL_SUBTEST_7( qr_invertible<MatrixXcf>() );
 | 
|---|
| 115 |     CALL_SUBTEST_8( qr_invertible<MatrixXcd>() );
 | 
|---|
| 116 |   }
 | 
|---|
| 117 | 
 | 
|---|
| 118 |   CALL_SUBTEST_9(qr_verify_assert<Matrix3f>());
 | 
|---|
| 119 |   CALL_SUBTEST_10(qr_verify_assert<Matrix3d>());
 | 
|---|
| 120 |   CALL_SUBTEST_1(qr_verify_assert<MatrixXf>());
 | 
|---|
| 121 |   CALL_SUBTEST_6(qr_verify_assert<MatrixXd>());
 | 
|---|
| 122 |   CALL_SUBTEST_7(qr_verify_assert<MatrixXcf>());
 | 
|---|
| 123 |   CALL_SUBTEST_8(qr_verify_assert<MatrixXcd>());
 | 
|---|
| 124 | 
 | 
|---|
| 125 |   // Test problem size constructors
 | 
|---|
| 126 |   CALL_SUBTEST_12(HouseholderQR<MatrixXf>(10, 20));
 | 
|---|
| 127 | }
 | 
|---|