aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/qr_fullpivoting.cpp
blob: 1e0601a83451d868866f624e323b0bc4074e68d3 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 3 of the License, or (at your option) any later version.
//
// Alternatively, you can redistribute it and/or
// modify it under the terms of the GNU General Public License as
// published by the Free Software Foundation; either version 2 of
// the License, or (at your option) any later version.
//
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.

#include "main.h"
#include <Eigen/QR>

template<typename MatrixType> void qr()
{
  /* this test covers the following files: QR.h */
  int rows = ei_random<int>(20,200), cols = ei_random<int>(20,200), cols2 = ei_random<int>(20,200);
  int rank = ei_random<int>(1, std::min(rows, cols)-1);

  typedef typename MatrixType::Scalar Scalar;
  typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> SquareMatrixType;
  typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
  MatrixType m1;
  createRandomMatrixOfRank(rank,rows,cols,m1);
  FullPivotingHouseholderQR<MatrixType> qr(m1);
  VERIFY_IS_APPROX(rank, qr.rank());
  VERIFY(cols - qr.rank() == qr.dimensionOfKernel());
  VERIFY(!qr.isInjective());
  VERIFY(!qr.isInvertible());
  VERIFY(!qr.isSurjective());

  
  MatrixType r = qr.matrixQR();
  // FIXME need better way to construct trapezoid
  for(int i = 0; i < rows; i++) for(int j = 0; j < cols; j++) if(i>j) r(i,j) = Scalar(0);
  
  MatrixType b = qr.matrixQ() * r;

  MatrixType c = MatrixType::Zero(rows,cols);
  
  for(int i = 0; i < cols; ++i) c.col(qr.colsPermutation().coeff(i)) = b.col(i);
  VERIFY_IS_APPROX(m1, c);
  
  MatrixType m2 = MatrixType::Random(cols,cols2);
  MatrixType m3 = m1*m2;
  m2 = MatrixType::Random(cols,cols2);
  VERIFY(qr.solve(m3, &m2));
  VERIFY_IS_APPROX(m3, m1*m2);
  m3 = MatrixType::Random(rows,cols2);
  VERIFY(!qr.solve(m3, &m2));
}

template<typename MatrixType> void qr_invertible()
{
  typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
  int size = ei_random<int>(10,50);

  MatrixType m1(size, size), m2(size, size), m3(size, size);
  m1 = MatrixType::Random(size,size);

  if (ei_is_same_type<RealScalar,float>::ret)
  {
    // let's build a matrix more stable to inverse
    MatrixType a = MatrixType::Random(size,size*2);
    m1 += a * a.adjoint();
  }

  FullPivotingHouseholderQR<MatrixType> qr(m1);
  VERIFY(qr.isInjective());
  VERIFY(qr.isInvertible());
  VERIFY(qr.isSurjective());

  m3 = MatrixType::Random(size,size);
  VERIFY(qr.solve(m3, &m2));
  VERIFY_IS_APPROX(m3, m1*m2);
}

template<typename MatrixType> void qr_verify_assert()
{
  MatrixType tmp;

  FullPivotingHouseholderQR<MatrixType> qr;
  VERIFY_RAISES_ASSERT(qr.matrixR())
  VERIFY_RAISES_ASSERT(qr.solve(tmp,&tmp))
  VERIFY_RAISES_ASSERT(qr.matrixQ())
}

void test_qr_fullpivoting()
{
 for(int i = 0; i < 1; i++) {
    // FIXME : very weird bug here
//     CALL_SUBTEST( qr(Matrix2f()) );
    CALL_SUBTEST( qr<MatrixXf>() );
    CALL_SUBTEST( qr<MatrixXd>() );
    CALL_SUBTEST( qr<MatrixXcd>() );
  }

  for(int i = 0; i < g_repeat; i++) {
    CALL_SUBTEST( qr_invertible<MatrixXf>() );
    CALL_SUBTEST( qr_invertible<MatrixXd>() );
    CALL_SUBTEST( qr_invertible<MatrixXcf>() );
    CALL_SUBTEST( qr_invertible<MatrixXcd>() );
  }

  CALL_SUBTEST(qr_verify_assert<Matrix3f>());
  CALL_SUBTEST(qr_verify_assert<Matrix3d>());
  CALL_SUBTEST(qr_verify_assert<MatrixXf>());
  CALL_SUBTEST(qr_verify_assert<MatrixXd>());
  CALL_SUBTEST(qr_verify_assert<MatrixXcf>());
  CALL_SUBTEST(qr_verify_assert<MatrixXcd>());
}