qr_fullpivoting.cpp
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
00006 //
00007 // Eigen is free software; you can redistribute it and/or
00008 // modify it under the terms of the GNU Lesser General Public
00009 // License as published by the Free Software Foundation; either
00010 // version 3 of the License, or (at your option) any later version.
00011 //
00012 // Alternatively, you can redistribute it and/or
00013 // modify it under the terms of the GNU General Public License as
00014 // published by the Free Software Foundation; either version 2 of
00015 // the License, or (at your option) any later version.
00016 //
00017 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00018 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00019 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00020 // GNU General Public License for more details.
00021 //
00022 // You should have received a copy of the GNU Lesser General Public
00023 // License and a copy of the GNU General Public License along with
00024 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00025 
00026 #include "main.h"
00027 #include <Eigen/QR>
00028 
00029 template<typename MatrixType> void qr()
00030 {
00031   typedef typename MatrixType::Index Index;
00032 
00033   Index rows = internal::random<Index>(20,200), cols = internal::random<int>(20,200), cols2 = internal::random<int>(20,200);
00034   Index rank = internal::random<Index>(1, std::min(rows, cols)-1);
00035 
00036   typedef typename MatrixType::Scalar Scalar;
00037   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
00038   typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
00039   MatrixType m1;
00040   createRandomPIMatrixOfRank(rank,rows,cols,m1);
00041   FullPivHouseholderQR<MatrixType> qr(m1);
00042   VERIFY(rank == qr.rank());
00043   VERIFY(cols - qr.rank() == qr.dimensionOfKernel());
00044   VERIFY(!qr.isInjective());
00045   VERIFY(!qr.isInvertible());
00046   VERIFY(!qr.isSurjective());
00047 
00048   MatrixType r = qr.matrixQR();
00049   
00050   MatrixQType q = qr.matrixQ();
00051   VERIFY_IS_UNITARY(q);
00052   
00053   // FIXME need better way to construct trapezoid
00054   for(int i = 0; i < rows; i++) for(int j = 0; j < cols; j++) if(i>j) r(i,j) = Scalar(0);
00055 
00056   MatrixType c = qr.matrixQ() * r * qr.colsPermutation().inverse();
00057 
00058   VERIFY_IS_APPROX(m1, c);
00059 
00060   MatrixType m2 = MatrixType::Random(cols,cols2);
00061   MatrixType m3 = m1*m2;
00062   m2 = MatrixType::Random(cols,cols2);
00063   m2 = qr.solve(m3);
00064   VERIFY_IS_APPROX(m3, m1*m2);
00065 }
00066 
00067 template<typename MatrixType> void qr_invertible()
00068 {
00069   typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
00070   typedef typename MatrixType::Scalar Scalar;
00071 
00072   int size = internal::random<int>(10,50);
00073 
00074   MatrixType m1(size, size), m2(size, size), m3(size, size);
00075   m1 = MatrixType::Random(size,size);
00076 
00077   if (internal::is_same<RealScalar,float>::value)
00078   {
00079     // let's build a matrix more stable to inverse
00080     MatrixType a = MatrixType::Random(size,size*2);
00081     m1 += a * a.adjoint();
00082   }
00083 
00084   FullPivHouseholderQR<MatrixType> qr(m1);
00085   VERIFY(qr.isInjective());
00086   VERIFY(qr.isInvertible());
00087   VERIFY(qr.isSurjective());
00088 
00089   m3 = MatrixType::Random(size,size);
00090   m2 = qr.solve(m3);
00091   VERIFY_IS_APPROX(m3, m1*m2);
00092 
00093   // now construct a matrix with prescribed determinant
00094   m1.setZero();
00095   for(int i = 0; i < size; i++) m1(i,i) = internal::random<Scalar>();
00096   RealScalar absdet = internal::abs(m1.diagonal().prod());
00097   m3 = qr.matrixQ(); // get a unitary
00098   m1 = m3 * m1 * m3;
00099   qr.compute(m1);
00100   VERIFY_IS_APPROX(absdet, qr.absDeterminant());
00101   VERIFY_IS_APPROX(internal::log(absdet), qr.logAbsDeterminant());
00102 }
00103 
00104 template<typename MatrixType> void qr_verify_assert()
00105 {
00106   MatrixType tmp;
00107 
00108   FullPivHouseholderQR<MatrixType> qr;
00109   VERIFY_RAISES_ASSERT(qr.matrixQR())
00110   VERIFY_RAISES_ASSERT(qr.solve(tmp))
00111   VERIFY_RAISES_ASSERT(qr.matrixQ())
00112   VERIFY_RAISES_ASSERT(qr.dimensionOfKernel())
00113   VERIFY_RAISES_ASSERT(qr.isInjective())
00114   VERIFY_RAISES_ASSERT(qr.isSurjective())
00115   VERIFY_RAISES_ASSERT(qr.isInvertible())
00116   VERIFY_RAISES_ASSERT(qr.inverse())
00117   VERIFY_RAISES_ASSERT(qr.absDeterminant())
00118   VERIFY_RAISES_ASSERT(qr.logAbsDeterminant())
00119 }
00120 
00121 void test_qr_fullpivoting()
00122 {
00123  for(int i = 0; i < 1; i++) {
00124     // FIXME : very weird bug here
00125 //     CALL_SUBTEST(qr(Matrix2f()) );
00126     CALL_SUBTEST_1( qr<MatrixXf>() );
00127     CALL_SUBTEST_2( qr<MatrixXd>() );
00128     CALL_SUBTEST_3( qr<MatrixXcd>() );
00129   }
00130 
00131   for(int i = 0; i < g_repeat; i++) {
00132     CALL_SUBTEST_1( qr_invertible<MatrixXf>() );
00133     CALL_SUBTEST_2( qr_invertible<MatrixXd>() );
00134     CALL_SUBTEST_4( qr_invertible<MatrixXcf>() );
00135     CALL_SUBTEST_3( qr_invertible<MatrixXcd>() );
00136   }
00137 
00138   CALL_SUBTEST_5(qr_verify_assert<Matrix3f>());
00139   CALL_SUBTEST_6(qr_verify_assert<Matrix3d>());
00140   CALL_SUBTEST_1(qr_verify_assert<MatrixXf>());
00141   CALL_SUBTEST_2(qr_verify_assert<MatrixXd>());
00142   CALL_SUBTEST_4(qr_verify_assert<MatrixXcf>());
00143   CALL_SUBTEST_3(qr_verify_assert<MatrixXcd>());
00144 
00145   // Test problem size constructors
00146   CALL_SUBTEST_7(FullPivHouseholderQR<MatrixXf>(10, 20));
00147 }


re_vision
Author(s): Dorian Galvez-Lopez
autogenerated on Sun Jan 5 2014 11:32:16