2009-08-22 13:13:21 +08:00
|
|
|
// This file is part of Eigen, a lightweight C++ template library
|
|
|
|
// for linear algebra.
|
|
|
|
//
|
2010-06-25 05:21:58 +08:00
|
|
|
// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
|
2009-08-22 13:13:21 +08:00
|
|
|
// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
|
|
|
|
//
|
2012-07-14 02:42:47 +08:00
|
|
|
// This Source Code Form is subject to the terms of the Mozilla
|
|
|
|
// Public License v. 2.0. If a copy of the MPL was not distributed
|
|
|
|
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
2009-08-22 13:13:21 +08:00
|
|
|
|
|
|
|
#include "main.h"
|
|
|
|
#include <Eigen/QR>
|
|
|
|
|
|
|
|
template<typename MatrixType> void qr()
|
|
|
|
{
|
2010-06-20 23:37:56 +08:00
|
|
|
typedef typename MatrixType::Index Index;
|
|
|
|
|
2016-07-20 21:14:20 +08:00
|
|
|
Index max_size = EIGEN_TEST_MAX_SIZE;
|
|
|
|
Index min_size = numext::maxi(1,EIGEN_TEST_MAX_SIZE/10);
|
|
|
|
Index rows = internal::random<Index>(min_size,max_size),
|
|
|
|
cols = internal::random<Index>(min_size,max_size),
|
|
|
|
cols2 = internal::random<Index>(min_size,max_size),
|
|
|
|
rank = internal::random<Index>(1, (std::min)(rows, cols)-1);
|
2009-08-22 13:13:21 +08:00
|
|
|
|
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
2009-10-13 10:33:51 +08:00
|
|
|
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
|
2009-08-22 13:13:21 +08:00
|
|
|
MatrixType m1;
|
2010-02-24 04:40:24 +08:00
|
|
|
createRandomPIMatrixOfRank(rank,rows,cols,m1);
|
2009-10-29 06:19:29 +08:00
|
|
|
FullPivHouseholderQR<MatrixType> qr(m1);
|
2015-06-26 22:22:49 +08:00
|
|
|
VERIFY_IS_EQUAL(rank, qr.rank());
|
|
|
|
VERIFY_IS_EQUAL(cols - qr.rank(), qr.dimensionOfKernel());
|
2009-08-24 12:23:35 +08:00
|
|
|
VERIFY(!qr.isInjective());
|
|
|
|
VERIFY(!qr.isInvertible());
|
|
|
|
VERIFY(!qr.isSurjective());
|
|
|
|
|
2009-08-22 13:13:21 +08:00
|
|
|
MatrixType r = qr.matrixQR();
|
2009-10-13 10:33:51 +08:00
|
|
|
|
|
|
|
MatrixQType q = qr.matrixQ();
|
|
|
|
VERIFY_IS_UNITARY(q);
|
|
|
|
|
2009-08-22 13:13:21 +08:00
|
|
|
// 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);
|
2009-09-16 21:56:20 +08:00
|
|
|
|
2009-11-17 21:14:54 +08:00
|
|
|
MatrixType c = qr.matrixQ() * r * qr.colsPermutation().inverse();
|
2009-08-22 13:13:21 +08:00
|
|
|
|
|
|
|
VERIFY_IS_APPROX(m1, c);
|
2014-07-21 17:45:54 +08:00
|
|
|
|
|
|
|
// stress the ReturnByValue mechanism
|
|
|
|
MatrixType tmp;
|
|
|
|
VERIFY_IS_APPROX(tmp.noalias() = qr.matrixQ() * r, (qr.matrixQ() * r).eval());
|
|
|
|
|
2009-08-22 13:13:21 +08:00
|
|
|
MatrixType m2 = MatrixType::Random(cols,cols2);
|
|
|
|
MatrixType m3 = m1*m2;
|
|
|
|
m2 = MatrixType::Random(cols,cols2);
|
2009-11-08 23:21:26 +08:00
|
|
|
m2 = qr.solve(m3);
|
2009-08-22 13:13:21 +08:00
|
|
|
VERIFY_IS_APPROX(m3, m1*m2);
|
2016-10-06 15:55:50 +08:00
|
|
|
|
|
|
|
{
|
|
|
|
Index size = rows;
|
|
|
|
do {
|
|
|
|
m1 = MatrixType::Random(size,size);
|
|
|
|
qr.compute(m1);
|
|
|
|
} while(!qr.isInvertible());
|
|
|
|
MatrixType m1_inv = qr.inverse();
|
|
|
|
m3 = m1 * MatrixType::Random(size,cols2);
|
|
|
|
m2 = qr.solve(m3);
|
|
|
|
VERIFY_IS_APPROX(m2, m1_inv*m3);
|
|
|
|
}
|
2009-08-22 13:13:21 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
template<typename MatrixType> void qr_invertible()
|
|
|
|
{
|
2012-11-06 22:25:50 +08:00
|
|
|
using std::log;
|
|
|
|
using std::abs;
|
2009-08-22 13:13:21 +08:00
|
|
|
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
|
2009-08-24 12:35:42 +08:00
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
|
|
|
|
2016-07-20 21:14:20 +08:00
|
|
|
Index max_size = numext::mini(50,EIGEN_TEST_MAX_SIZE);
|
|
|
|
Index min_size = numext::maxi(1,EIGEN_TEST_MAX_SIZE/10);
|
|
|
|
Index size = internal::random<Index>(min_size,max_size);
|
2009-08-22 13:13:21 +08:00
|
|
|
|
|
|
|
MatrixType m1(size, size), m2(size, size), m3(size, size);
|
|
|
|
m1 = MatrixType::Random(size,size);
|
|
|
|
|
2010-10-26 04:13:49 +08:00
|
|
|
if (internal::is_same<RealScalar,float>::value)
|
2009-08-22 13:13:21 +08:00
|
|
|
{
|
|
|
|
// let's build a matrix more stable to inverse
|
|
|
|
MatrixType a = MatrixType::Random(size,size*2);
|
|
|
|
m1 += a * a.adjoint();
|
|
|
|
}
|
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
FullPivHouseholderQR<MatrixType> qr(m1);
|
2009-08-24 12:23:35 +08:00
|
|
|
VERIFY(qr.isInjective());
|
|
|
|
VERIFY(qr.isInvertible());
|
|
|
|
VERIFY(qr.isSurjective());
|
|
|
|
|
2009-08-22 13:13:21 +08:00
|
|
|
m3 = MatrixType::Random(size,size);
|
2009-11-08 23:21:26 +08:00
|
|
|
m2 = qr.solve(m3);
|
2009-08-22 13:13:21 +08:00
|
|
|
VERIFY_IS_APPROX(m3, m1*m2);
|
2009-09-16 21:56:20 +08:00
|
|
|
|
2009-08-24 12:35:42 +08:00
|
|
|
// now construct a matrix with prescribed determinant
|
|
|
|
m1.setZero();
|
2010-10-25 22:15:22 +08:00
|
|
|
for(int i = 0; i < size; i++) m1(i,i) = internal::random<Scalar>();
|
2012-11-06 22:25:50 +08:00
|
|
|
RealScalar absdet = abs(m1.diagonal().prod());
|
2009-08-24 12:35:42 +08:00
|
|
|
m3 = qr.matrixQ(); // get a unitary
|
|
|
|
m1 = m3 * m1 * m3;
|
|
|
|
qr.compute(m1);
|
|
|
|
VERIFY_IS_APPROX(absdet, qr.absDeterminant());
|
2012-11-06 22:25:50 +08:00
|
|
|
VERIFY_IS_APPROX(log(absdet), qr.logAbsDeterminant());
|
2009-08-22 13:13:21 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
template<typename MatrixType> void qr_verify_assert()
|
|
|
|
{
|
|
|
|
MatrixType tmp;
|
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
FullPivHouseholderQR<MatrixType> qr;
|
2009-08-25 01:46:14 +08:00
|
|
|
VERIFY_RAISES_ASSERT(qr.matrixQR())
|
2009-11-08 23:21:26 +08:00
|
|
|
VERIFY_RAISES_ASSERT(qr.solve(tmp))
|
2009-08-22 13:13:21 +08:00
|
|
|
VERIFY_RAISES_ASSERT(qr.matrixQ())
|
2009-08-24 23:11:41 +08:00
|
|
|
VERIFY_RAISES_ASSERT(qr.dimensionOfKernel())
|
|
|
|
VERIFY_RAISES_ASSERT(qr.isInjective())
|
|
|
|
VERIFY_RAISES_ASSERT(qr.isSurjective())
|
|
|
|
VERIFY_RAISES_ASSERT(qr.isInvertible())
|
|
|
|
VERIFY_RAISES_ASSERT(qr.inverse())
|
2009-08-25 01:46:14 +08:00
|
|
|
VERIFY_RAISES_ASSERT(qr.absDeterminant())
|
|
|
|
VERIFY_RAISES_ASSERT(qr.logAbsDeterminant())
|
2009-08-22 13:13:21 +08:00
|
|
|
}
|
|
|
|
|
2009-08-24 06:04:33 +08:00
|
|
|
void test_qr_fullpivoting()
|
2009-08-22 13:13:21 +08:00
|
|
|
{
|
|
|
|
for(int i = 0; i < 1; i++) {
|
|
|
|
// FIXME : very weird bug here
|
2009-10-29 06:19:29 +08:00
|
|
|
// CALL_SUBTEST(qr(Matrix2f()) );
|
|
|
|
CALL_SUBTEST_1( qr<MatrixXf>() );
|
|
|
|
CALL_SUBTEST_2( qr<MatrixXd>() );
|
|
|
|
CALL_SUBTEST_3( qr<MatrixXcd>() );
|
2009-08-22 13:13:21 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
for(int i = 0; i < g_repeat; i++) {
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_1( qr_invertible<MatrixXf>() );
|
|
|
|
CALL_SUBTEST_2( qr_invertible<MatrixXd>() );
|
|
|
|
CALL_SUBTEST_4( qr_invertible<MatrixXcf>() );
|
|
|
|
CALL_SUBTEST_3( qr_invertible<MatrixXcd>() );
|
2009-08-22 13:13:21 +08:00
|
|
|
}
|
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_5(qr_verify_assert<Matrix3f>());
|
|
|
|
CALL_SUBTEST_6(qr_verify_assert<Matrix3d>());
|
|
|
|
CALL_SUBTEST_1(qr_verify_assert<MatrixXf>());
|
|
|
|
CALL_SUBTEST_2(qr_verify_assert<MatrixXd>());
|
|
|
|
CALL_SUBTEST_4(qr_verify_assert<MatrixXcf>());
|
|
|
|
CALL_SUBTEST_3(qr_verify_assert<MatrixXcd>());
|
2010-04-21 23:15:57 +08:00
|
|
|
|
|
|
|
// Test problem size constructors
|
|
|
|
CALL_SUBTEST_7(FullPivHouseholderQR<MatrixXf>(10, 20));
|
2013-11-19 19:53:46 +08:00
|
|
|
CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,10,20> >(10,20)));
|
|
|
|
CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,10,20> >(Matrix<float,10,20>::Random())));
|
|
|
|
CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,20,10> >(20,10)));
|
|
|
|
CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,20,10> >(Matrix<float,20,10>::Random())));
|
2009-08-22 13:13:21 +08:00
|
|
|
}
|