qr_colpivoting.cpp revision c981c48f5bc9aefeffc0bcb0cc3934c2fae179dd
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
14template<typename MatrixType> void qr()
15{
16  typedef typename MatrixType::Index Index;
17
18  Index rows = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE), cols = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE), cols2 = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE);
19  Index rank = internal::random<Index>(1, (std::min)(rows, cols)-1);
20
21  typedef typename MatrixType::Scalar Scalar;
22  typedef typename MatrixType::RealScalar RealScalar;
23  typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
24  typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
25  MatrixType m1;
26  createRandomPIMatrixOfRank(rank,rows,cols,m1);
27  ColPivHouseholderQR<MatrixType> qr(m1);
28  VERIFY(rank == qr.rank());
29  VERIFY(cols - qr.rank() == qr.dimensionOfKernel());
30  VERIFY(!qr.isInjective());
31  VERIFY(!qr.isInvertible());
32  VERIFY(!qr.isSurjective());
33
34  MatrixQType q = qr.householderQ();
35  VERIFY_IS_UNITARY(q);
36
37  MatrixType r = qr.matrixQR().template triangularView<Upper>();
38  MatrixType c = q * r * qr.colsPermutation().inverse();
39  VERIFY_IS_APPROX(m1, c);
40
41  MatrixType m2 = MatrixType::Random(cols,cols2);
42  MatrixType m3 = m1*m2;
43  m2 = MatrixType::Random(cols,cols2);
44  m2 = qr.solve(m3);
45  VERIFY_IS_APPROX(m3, m1*m2);
46}
47
48template<typename MatrixType, int Cols2> void qr_fixedsize()
49{
50  enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime };
51  typedef typename MatrixType::Scalar Scalar;
52  int rank = internal::random<int>(1, (std::min)(int(Rows), int(Cols))-1);
53  Matrix<Scalar,Rows,Cols> m1;
54  createRandomPIMatrixOfRank(rank,Rows,Cols,m1);
55  ColPivHouseholderQR<Matrix<Scalar,Rows,Cols> > qr(m1);
56  VERIFY(rank == qr.rank());
57  VERIFY(Cols - qr.rank() == qr.dimensionOfKernel());
58  VERIFY(qr.isInjective() == (rank == Rows));
59  VERIFY(qr.isSurjective() == (rank == Cols));
60  VERIFY(qr.isInvertible() == (qr.isInjective() && qr.isSurjective()));
61
62  Matrix<Scalar,Rows,Cols> r = qr.matrixQR().template triangularView<Upper>();
63  Matrix<Scalar,Rows,Cols> c = qr.householderQ() * r * qr.colsPermutation().inverse();
64  VERIFY_IS_APPROX(m1, c);
65
66  Matrix<Scalar,Cols,Cols2> m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2);
67  Matrix<Scalar,Rows,Cols2> m3 = m1*m2;
68  m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2);
69  m2 = qr.solve(m3);
70  VERIFY_IS_APPROX(m3, m1*m2);
71}
72
73template<typename MatrixType> void qr_invertible()
74{
75  typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
76  typedef typename MatrixType::Scalar Scalar;
77
78  int size = internal::random<int>(10,50);
79
80  MatrixType m1(size, size), m2(size, size), m3(size, size);
81  m1 = MatrixType::Random(size,size);
82
83  if (internal::is_same<RealScalar,float>::value)
84  {
85    // let's build a matrix more stable to inverse
86    MatrixType a = MatrixType::Random(size,size*2);
87    m1 += a * a.adjoint();
88  }
89
90  ColPivHouseholderQR<MatrixType> qr(m1);
91  m3 = MatrixType::Random(size,size);
92  m2 = qr.solve(m3);
93  //VERIFY_IS_APPROX(m3, m1*m2);
94
95  // now construct a matrix with prescribed determinant
96  m1.setZero();
97  for(int i = 0; i < size; i++) m1(i,i) = internal::random<Scalar>();
98  RealScalar absdet = internal::abs(m1.diagonal().prod());
99  m3 = qr.householderQ(); // get a unitary
100  m1 = m3 * m1 * m3;
101  qr.compute(m1);
102  VERIFY_IS_APPROX(absdet, qr.absDeterminant());
103  VERIFY_IS_APPROX(internal::log(absdet), qr.logAbsDeterminant());
104}
105
106template<typename MatrixType> void qr_verify_assert()
107{
108  MatrixType tmp;
109
110  ColPivHouseholderQR<MatrixType> qr;
111  VERIFY_RAISES_ASSERT(qr.matrixQR())
112  VERIFY_RAISES_ASSERT(qr.solve(tmp))
113  VERIFY_RAISES_ASSERT(qr.householderQ())
114  VERIFY_RAISES_ASSERT(qr.dimensionOfKernel())
115  VERIFY_RAISES_ASSERT(qr.isInjective())
116  VERIFY_RAISES_ASSERT(qr.isSurjective())
117  VERIFY_RAISES_ASSERT(qr.isInvertible())
118  VERIFY_RAISES_ASSERT(qr.inverse())
119  VERIFY_RAISES_ASSERT(qr.absDeterminant())
120  VERIFY_RAISES_ASSERT(qr.logAbsDeterminant())
121}
122
123void test_qr_colpivoting()
124{
125  for(int i = 0; i < g_repeat; i++) {
126    CALL_SUBTEST_1( qr<MatrixXf>() );
127    CALL_SUBTEST_2( qr<MatrixXd>() );
128    CALL_SUBTEST_3( qr<MatrixXcd>() );
129    CALL_SUBTEST_4(( qr_fixedsize<Matrix<float,3,5>, 4 >() ));
130    CALL_SUBTEST_5(( qr_fixedsize<Matrix<double,6,2>, 3 >() ));
131    CALL_SUBTEST_5(( qr_fixedsize<Matrix<double,1,1>, 1 >() ));
132  }
133
134  for(int i = 0; i < g_repeat; i++) {
135    CALL_SUBTEST_1( qr_invertible<MatrixXf>() );
136    CALL_SUBTEST_2( qr_invertible<MatrixXd>() );
137    CALL_SUBTEST_6( qr_invertible<MatrixXcf>() );
138    CALL_SUBTEST_3( qr_invertible<MatrixXcd>() );
139  }
140
141  CALL_SUBTEST_7(qr_verify_assert<Matrix3f>());
142  CALL_SUBTEST_8(qr_verify_assert<Matrix3d>());
143  CALL_SUBTEST_1(qr_verify_assert<MatrixXf>());
144  CALL_SUBTEST_2(qr_verify_assert<MatrixXd>());
145  CALL_SUBTEST_6(qr_verify_assert<MatrixXcf>());
146  CALL_SUBTEST_3(qr_verify_assert<MatrixXcd>());
147
148  // Test problem size constructors
149  CALL_SUBTEST_9(ColPivHouseholderQR<MatrixXf>(10, 20));
150}
151