• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
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 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #include "main.h"
12 #include <Eigen/QR>
13 
qr()14 template<typename MatrixType> void qr()
15 {
16   typedef typename MatrixType::Index Index;
17 
18   Index max_size = EIGEN_TEST_MAX_SIZE;
19   Index min_size = numext::maxi(1,EIGEN_TEST_MAX_SIZE/10);
20   Index rows  = internal::random<Index>(min_size,max_size),
21         cols  = internal::random<Index>(min_size,max_size),
22         cols2 = internal::random<Index>(min_size,max_size),
23         rank  = internal::random<Index>(1, (std::min)(rows, cols)-1);
24 
25   typedef typename MatrixType::Scalar Scalar;
26   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
27   MatrixType m1;
28   createRandomPIMatrixOfRank(rank,rows,cols,m1);
29   FullPivHouseholderQR<MatrixType> qr(m1);
30   VERIFY_IS_EQUAL(rank, qr.rank());
31   VERIFY_IS_EQUAL(cols - qr.rank(), qr.dimensionOfKernel());
32   VERIFY(!qr.isInjective());
33   VERIFY(!qr.isInvertible());
34   VERIFY(!qr.isSurjective());
35 
36   MatrixType r = qr.matrixQR();
37 
38   MatrixQType q = qr.matrixQ();
39   VERIFY_IS_UNITARY(q);
40 
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   MatrixType c = qr.matrixQ() * r * qr.colsPermutation().inverse();
45 
46   VERIFY_IS_APPROX(m1, c);
47 
48   // stress the ReturnByValue mechanism
49   MatrixType tmp;
50   VERIFY_IS_APPROX(tmp.noalias() = qr.matrixQ() * r, (qr.matrixQ() * r).eval());
51 
52   MatrixType m2 = MatrixType::Random(cols,cols2);
53   MatrixType m3 = m1*m2;
54   m2 = MatrixType::Random(cols,cols2);
55   m2 = qr.solve(m3);
56   VERIFY_IS_APPROX(m3, m1*m2);
57 
58   {
59     Index size = rows;
60     do {
61       m1 = MatrixType::Random(size,size);
62       qr.compute(m1);
63     } while(!qr.isInvertible());
64     MatrixType m1_inv = qr.inverse();
65     m3 = m1 * MatrixType::Random(size,cols2);
66     m2 = qr.solve(m3);
67     VERIFY_IS_APPROX(m2, m1_inv*m3);
68   }
69 }
70 
qr_invertible()71 template<typename MatrixType> void qr_invertible()
72 {
73   using std::log;
74   using std::abs;
75   typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
76   typedef typename MatrixType::Scalar Scalar;
77 
78   Index max_size = numext::mini(50,EIGEN_TEST_MAX_SIZE);
79   Index min_size = numext::maxi(1,EIGEN_TEST_MAX_SIZE/10);
80   Index size = internal::random<Index>(min_size,max_size);
81 
82   MatrixType m1(size, size), m2(size, size), m3(size, size);
83   m1 = MatrixType::Random(size,size);
84 
85   if (internal::is_same<RealScalar,float>::value)
86   {
87     // let's build a matrix more stable to inverse
88     MatrixType a = MatrixType::Random(size,size*2);
89     m1 += a * a.adjoint();
90   }
91 
92   FullPivHouseholderQR<MatrixType> qr(m1);
93   VERIFY(qr.isInjective());
94   VERIFY(qr.isInvertible());
95   VERIFY(qr.isSurjective());
96 
97   m3 = MatrixType::Random(size,size);
98   m2 = qr.solve(m3);
99   VERIFY_IS_APPROX(m3, m1*m2);
100 
101   // now construct a matrix with prescribed determinant
102   m1.setZero();
103   for(int i = 0; i < size; i++) m1(i,i) = internal::random<Scalar>();
104   RealScalar absdet = abs(m1.diagonal().prod());
105   m3 = qr.matrixQ(); // get a unitary
106   m1 = m3 * m1 * m3;
107   qr.compute(m1);
108   VERIFY_IS_APPROX(absdet, qr.absDeterminant());
109   VERIFY_IS_APPROX(log(absdet), qr.logAbsDeterminant());
110 }
111 
qr_verify_assert()112 template<typename MatrixType> void qr_verify_assert()
113 {
114   MatrixType tmp;
115 
116   FullPivHouseholderQR<MatrixType> qr;
117   VERIFY_RAISES_ASSERT(qr.matrixQR())
118   VERIFY_RAISES_ASSERT(qr.solve(tmp))
119   VERIFY_RAISES_ASSERT(qr.matrixQ())
120   VERIFY_RAISES_ASSERT(qr.dimensionOfKernel())
121   VERIFY_RAISES_ASSERT(qr.isInjective())
122   VERIFY_RAISES_ASSERT(qr.isSurjective())
123   VERIFY_RAISES_ASSERT(qr.isInvertible())
124   VERIFY_RAISES_ASSERT(qr.inverse())
125   VERIFY_RAISES_ASSERT(qr.absDeterminant())
126   VERIFY_RAISES_ASSERT(qr.logAbsDeterminant())
127 }
128 
test_qr_fullpivoting()129 void test_qr_fullpivoting()
130 {
131  for(int i = 0; i < 1; i++) {
132     // FIXME : very weird bug here
133 //     CALL_SUBTEST(qr(Matrix2f()) );
134     CALL_SUBTEST_1( qr<MatrixXf>() );
135     CALL_SUBTEST_2( qr<MatrixXd>() );
136     CALL_SUBTEST_3( qr<MatrixXcd>() );
137   }
138 
139   for(int i = 0; i < g_repeat; i++) {
140     CALL_SUBTEST_1( qr_invertible<MatrixXf>() );
141     CALL_SUBTEST_2( qr_invertible<MatrixXd>() );
142     CALL_SUBTEST_4( qr_invertible<MatrixXcf>() );
143     CALL_SUBTEST_3( qr_invertible<MatrixXcd>() );
144   }
145 
146   CALL_SUBTEST_5(qr_verify_assert<Matrix3f>());
147   CALL_SUBTEST_6(qr_verify_assert<Matrix3d>());
148   CALL_SUBTEST_1(qr_verify_assert<MatrixXf>());
149   CALL_SUBTEST_2(qr_verify_assert<MatrixXd>());
150   CALL_SUBTEST_4(qr_verify_assert<MatrixXcf>());
151   CALL_SUBTEST_3(qr_verify_assert<MatrixXcd>());
152 
153   // Test problem size constructors
154   CALL_SUBTEST_7(FullPivHouseholderQR<MatrixXf>(10, 20));
155   CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,10,20> >(10,20)));
156   CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,10,20> >(Matrix<float,10,20>::Random())));
157   CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,20,10> >(20,10)));
158   CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,20,10> >(Matrix<float,20,10>::Random())));
159 }
160