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 13template<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 typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType; 23 24 MatrixType a = MatrixType::Random(rows,cols); 25 HouseholderQR<MatrixType> qrOfA(a); 26 27 MatrixQType q = qrOfA.householderQ(); 28 VERIFY_IS_UNITARY(q); 29 30 MatrixType r = qrOfA.matrixQR().template triangularView<Upper>(); 31 VERIFY_IS_APPROX(a, qrOfA.householderQ() * r); 32} 33 34template<typename MatrixType, int Cols2> void qr_fixedsize() 35{ 36 enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime }; 37 typedef typename MatrixType::Scalar Scalar; 38 Matrix<Scalar,Rows,Cols> m1 = Matrix<Scalar,Rows,Cols>::Random(); 39 HouseholderQR<Matrix<Scalar,Rows,Cols> > qr(m1); 40 41 Matrix<Scalar,Rows,Cols> r = qr.matrixQR(); 42 // FIXME need better way to construct trapezoid 43 for(int i = 0; i < Rows; i++) for(int j = 0; j < Cols; j++) if(i>j) r(i,j) = Scalar(0); 44 45 VERIFY_IS_APPROX(m1, qr.householderQ() * r); 46 47 Matrix<Scalar,Cols,Cols2> m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2); 48 Matrix<Scalar,Rows,Cols2> m3 = m1*m2; 49 m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2); 50 m2 = qr.solve(m3); 51 VERIFY_IS_APPROX(m3, m1*m2); 52} 53 54template<typename MatrixType> void qr_invertible() 55{ 56 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar; 57 typedef typename MatrixType::Scalar Scalar; 58 59 int size = internal::random<int>(10,50); 60 61 MatrixType m1(size, size), m2(size, size), m3(size, size); 62 m1 = MatrixType::Random(size,size); 63 64 if (internal::is_same<RealScalar,float>::value) 65 { 66 // let's build a matrix more stable to inverse 67 MatrixType a = MatrixType::Random(size,size*2); 68 m1 += a * a.adjoint(); 69 } 70 71 HouseholderQR<MatrixType> qr(m1); 72 m3 = MatrixType::Random(size,size); 73 m2 = qr.solve(m3); 74 VERIFY_IS_APPROX(m3, m1*m2); 75 76 // now construct a matrix with prescribed determinant 77 m1.setZero(); 78 for(int i = 0; i < size; i++) m1(i,i) = internal::random<Scalar>(); 79 RealScalar absdet = internal::abs(m1.diagonal().prod()); 80 m3 = qr.householderQ(); // get a unitary 81 m1 = m3 * m1 * m3; 82 qr.compute(m1); 83 VERIFY_IS_APPROX(absdet, qr.absDeterminant()); 84 VERIFY_IS_APPROX(internal::log(absdet), qr.logAbsDeterminant()); 85} 86 87template<typename MatrixType> void qr_verify_assert() 88{ 89 MatrixType tmp; 90 91 HouseholderQR<MatrixType> qr; 92 VERIFY_RAISES_ASSERT(qr.matrixQR()) 93 VERIFY_RAISES_ASSERT(qr.solve(tmp)) 94 VERIFY_RAISES_ASSERT(qr.householderQ()) 95 VERIFY_RAISES_ASSERT(qr.absDeterminant()) 96 VERIFY_RAISES_ASSERT(qr.logAbsDeterminant()) 97} 98 99void test_qr() 100{ 101 for(int i = 0; i < g_repeat; i++) { 102 CALL_SUBTEST_1( qr(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); 103 CALL_SUBTEST_2( qr(MatrixXcd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2),internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))) ); 104 CALL_SUBTEST_3(( qr_fixedsize<Matrix<float,3,4>, 2 >() )); 105 CALL_SUBTEST_4(( qr_fixedsize<Matrix<double,6,2>, 4 >() )); 106 CALL_SUBTEST_5(( qr_fixedsize<Matrix<double,2,5>, 7 >() )); 107 CALL_SUBTEST_11( qr(Matrix<float,1,1>()) ); 108 } 109 110 for(int i = 0; i < g_repeat; i++) { 111 CALL_SUBTEST_1( qr_invertible<MatrixXf>() ); 112 CALL_SUBTEST_6( qr_invertible<MatrixXd>() ); 113 CALL_SUBTEST_7( qr_invertible<MatrixXcf>() ); 114 CALL_SUBTEST_8( qr_invertible<MatrixXcd>() ); 115 } 116 117 CALL_SUBTEST_9(qr_verify_assert<Matrix3f>()); 118 CALL_SUBTEST_10(qr_verify_assert<Matrix3d>()); 119 CALL_SUBTEST_1(qr_verify_assert<MatrixXf>()); 120 CALL_SUBTEST_6(qr_verify_assert<MatrixXd>()); 121 CALL_SUBTEST_7(qr_verify_assert<MatrixXcf>()); 122 CALL_SUBTEST_8(qr_verify_assert<MatrixXcd>()); 123 124 // Test problem size constructors 125 CALL_SUBTEST_12(HouseholderQR<MatrixXf>(10, 20)); 126} 127