Home | History | Annotate | Download | only in test
      1 // This file is part of Eigen, a lightweight C++ template library
      2 // for linear algebra.
      3 //
      4 // Copyright (C) 2009-2010 Benoit Jacob <jacob.benoit.1 (at) gmail.com>
      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 householder(const MatrixType& m)
     14 {
     15   typedef typename MatrixType::Index Index;
     16   static bool even = true;
     17   even = !even;
     18   /* this test covers the following files:
     19      Householder.h
     20   */
     21   Index rows = m.rows();
     22   Index cols = m.cols();
     23 
     24   typedef typename MatrixType::Scalar Scalar;
     25   typedef typename NumTraits<Scalar>::Real RealScalar;
     26   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
     27   typedef Matrix<Scalar, internal::decrement_size<MatrixType::RowsAtCompileTime>::ret, 1> EssentialVectorType;
     28   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
     29   typedef Matrix<Scalar, Dynamic, MatrixType::ColsAtCompileTime> HBlockMatrixType;
     30   typedef Matrix<Scalar, Dynamic, 1> HCoeffsVectorType;
     31 
     32   typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> RightSquareMatrixType;
     33   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, Dynamic> VBlockMatrixType;
     34   typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::RowsAtCompileTime> TMatrixType;
     35 
     36   Matrix<Scalar, EIGEN_SIZE_MAX(MatrixType::RowsAtCompileTime,MatrixType::ColsAtCompileTime), 1> _tmp((std::max)(rows,cols));
     37   Scalar* tmp = &_tmp.coeffRef(0,0);
     38 
     39   Scalar beta;
     40   RealScalar alpha;
     41   EssentialVectorType essential;
     42 
     43   VectorType v1 = VectorType::Random(rows), v2;
     44   v2 = v1;
     45   v1.makeHouseholder(essential, beta, alpha);
     46   v1.applyHouseholderOnTheLeft(essential,beta,tmp);
     47   VERIFY_IS_APPROX(v1.norm(), v2.norm());
     48   if(rows>=2) VERIFY_IS_MUCH_SMALLER_THAN(v1.tail(rows-1).norm(), v1.norm());
     49   v1 = VectorType::Random(rows);
     50   v2 = v1;
     51   v1.applyHouseholderOnTheLeft(essential,beta,tmp);
     52   VERIFY_IS_APPROX(v1.norm(), v2.norm());
     53 
     54   MatrixType m1(rows, cols),
     55              m2(rows, cols);
     56 
     57   v1 = VectorType::Random(rows);
     58   if(even) v1.tail(rows-1).setZero();
     59   m1.colwise() = v1;
     60   m2 = m1;
     61   m1.col(0).makeHouseholder(essential, beta, alpha);
     62   m1.applyHouseholderOnTheLeft(essential,beta,tmp);
     63   VERIFY_IS_APPROX(m1.norm(), m2.norm());
     64   if(rows>=2) VERIFY_IS_MUCH_SMALLER_THAN(m1.block(1,0,rows-1,cols).norm(), m1.norm());
     65   VERIFY_IS_MUCH_SMALLER_THAN(internal::imag(m1(0,0)), internal::real(m1(0,0)));
     66   VERIFY_IS_APPROX(internal::real(m1(0,0)), alpha);
     67 
     68   v1 = VectorType::Random(rows);
     69   if(even) v1.tail(rows-1).setZero();
     70   SquareMatrixType m3(rows,rows), m4(rows,rows);
     71   m3.rowwise() = v1.transpose();
     72   m4 = m3;
     73   m3.row(0).makeHouseholder(essential, beta, alpha);
     74   m3.applyHouseholderOnTheRight(essential,beta,tmp);
     75   VERIFY_IS_APPROX(m3.norm(), m4.norm());
     76   if(rows>=2) VERIFY_IS_MUCH_SMALLER_THAN(m3.block(0,1,rows,rows-1).norm(), m3.norm());
     77   VERIFY_IS_MUCH_SMALLER_THAN(internal::imag(m3(0,0)), internal::real(m3(0,0)));
     78   VERIFY_IS_APPROX(internal::real(m3(0,0)), alpha);
     79 
     80   // test householder sequence on the left with a shift
     81 
     82   Index shift = internal::random<Index>(0, std::max<Index>(rows-2,0));
     83   Index brows = rows - shift;
     84   m1.setRandom(rows, cols);
     85   HBlockMatrixType hbm = m1.block(shift,0,brows,cols);
     86   HouseholderQR<HBlockMatrixType> qr(hbm);
     87   m2 = m1;
     88   m2.block(shift,0,brows,cols) = qr.matrixQR();
     89   HCoeffsVectorType hc = qr.hCoeffs().conjugate();
     90   HouseholderSequence<MatrixType, HCoeffsVectorType> hseq(m2, hc);
     91   hseq.setLength(hc.size()).setShift(shift);
     92   VERIFY(hseq.length() == hc.size());
     93   VERIFY(hseq.shift() == shift);
     94 
     95   MatrixType m5 = m2;
     96   m5.block(shift,0,brows,cols).template triangularView<StrictlyLower>().setZero();
     97   VERIFY_IS_APPROX(hseq * m5, m1); // test applying hseq directly
     98   m3 = hseq;
     99   VERIFY_IS_APPROX(m3 * m5, m1); // test evaluating hseq to a dense matrix, then applying
    100 
    101   // test householder sequence on the right with a shift
    102 
    103   TMatrixType tm2 = m2.transpose();
    104   HouseholderSequence<TMatrixType, HCoeffsVectorType, OnTheRight> rhseq(tm2, hc);
    105   rhseq.setLength(hc.size()).setShift(shift);
    106   VERIFY_IS_APPROX(rhseq * m5, m1); // test applying rhseq directly
    107   m3 = rhseq;
    108   VERIFY_IS_APPROX(m3 * m5, m1); // test evaluating rhseq to a dense matrix, then applying
    109 }
    110 
    111 void test_householder()
    112 {
    113   for(int i = 0; i < g_repeat; i++) {
    114     CALL_SUBTEST_1( householder(Matrix<double,2,2>()) );
    115     CALL_SUBTEST_2( householder(Matrix<float,2,3>()) );
    116     CALL_SUBTEST_3( householder(Matrix<double,3,5>()) );
    117     CALL_SUBTEST_4( householder(Matrix<float,4,4>()) );
    118     CALL_SUBTEST_5( householder(MatrixXd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
    119     CALL_SUBTEST_6( householder(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
    120     CALL_SUBTEST_7( householder(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
    121     CALL_SUBTEST_8( householder(Matrix<double,1,1>()) );
    122   }
    123 }
    124