2009-08-22 13:13:21 +08:00
|
|
|
// 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()
|
|
|
|
{
|
2009-12-02 02:26:29 +08:00
|
|
|
int rows = ei_random<int>(2,200), cols = ei_random<int>(2,200), cols2 = ei_random<int>(2,200);
|
2009-08-22 13:13:21 +08:00
|
|
|
int rank = ei_random<int>(1, std::min(rows, cols)-1);
|
|
|
|
|
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
2009-12-02 02:26:29 +08:00
|
|
|
typedef typename MatrixType::RealScalar RealScalar;
|
2009-10-13 10:33:51 +08:00
|
|
|
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
|
2009-08-22 13:13:21 +08:00
|
|
|
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
|
|
|
|
MatrixType m1;
|
|
|
|
createRandomMatrixOfRank(rank,rows,cols,m1);
|
2009-10-29 06:19:29 +08:00
|
|
|
ColPivHouseholderQR<MatrixType> qr(m1);
|
2009-08-22 13:13:21 +08:00
|
|
|
VERIFY_IS_APPROX(rank, qr.rank());
|
2009-08-24 23:11:41 +08:00
|
|
|
VERIFY(cols - qr.rank() == qr.dimensionOfKernel());
|
|
|
|
VERIFY(!qr.isInjective());
|
|
|
|
VERIFY(!qr.isInvertible());
|
|
|
|
VERIFY(!qr.isSurjective());
|
2009-09-16 21:56:20 +08:00
|
|
|
|
2009-10-13 10:33:51 +08:00
|
|
|
MatrixQType q = qr.matrixQ();
|
|
|
|
VERIFY_IS_UNITARY(q);
|
2009-09-16 21:56:20 +08:00
|
|
|
|
2009-12-02 02:26:29 +08:00
|
|
|
MatrixType r = qr.matrixQR().template triangularView<UpperTriangular>();
|
|
|
|
MatrixType c = q * r * qr.colsPermutation().inverse();
|
2009-08-22 13:13:21 +08:00
|
|
|
VERIFY_IS_APPROX(m1, c);
|
2009-09-16 21:56:20 +08:00
|
|
|
|
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);
|
|
|
|
}
|
|
|
|
|
2009-09-28 21:40:18 +08:00
|
|
|
template<typename MatrixType, int Cols2> void qr_fixedsize()
|
|
|
|
{
|
|
|
|
enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime };
|
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
|
|
|
int rank = ei_random<int>(1, std::min(int(Rows), int(Cols))-1);
|
|
|
|
Matrix<Scalar,Rows,Cols> m1;
|
|
|
|
createRandomMatrixOfRank(rank,Rows,Cols,m1);
|
2009-10-29 06:19:29 +08:00
|
|
|
ColPivHouseholderQR<Matrix<Scalar,Rows,Cols> > qr(m1);
|
2009-09-28 21:40:18 +08:00
|
|
|
VERIFY_IS_APPROX(rank, qr.rank());
|
|
|
|
VERIFY(Cols - qr.rank() == qr.dimensionOfKernel());
|
|
|
|
VERIFY(!qr.isInjective());
|
|
|
|
VERIFY(!qr.isInvertible());
|
|
|
|
VERIFY(!qr.isSurjective());
|
|
|
|
|
2009-12-02 02:26:29 +08:00
|
|
|
Matrix<Scalar,Rows,Cols> r = qr.matrixQR().template triangularView<UpperTriangular>();
|
|
|
|
Matrix<Scalar,Rows,Cols> c = qr.matrixQ() * r * qr.colsPermutation().inverse();
|
2009-09-28 21:40:18 +08:00
|
|
|
VERIFY_IS_APPROX(m1, c);
|
|
|
|
|
|
|
|
Matrix<Scalar,Cols,Cols2> m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2);
|
|
|
|
Matrix<Scalar,Rows,Cols2> m3 = m1*m2;
|
|
|
|
m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2);
|
2009-11-08 23:21:26 +08:00
|
|
|
m2 = qr.solve(m3);
|
2009-09-28 21:40:18 +08:00
|
|
|
VERIFY_IS_APPROX(m3, m1*m2);
|
|
|
|
}
|
|
|
|
|
2009-08-22 13:13:21 +08:00
|
|
|
template<typename MatrixType> void qr_invertible()
|
|
|
|
{
|
|
|
|
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
|
2009-08-24 23:11:41 +08:00
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
|
|
|
|
2009-08-22 13:13:21 +08:00
|
|
|
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();
|
|
|
|
}
|
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
ColPivHouseholderQR<MatrixType> qr(m1);
|
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-12-02 02:26:29 +08:00
|
|
|
//VERIFY_IS_APPROX(m3, m1*m2);
|
2009-08-24 23:11:41 +08:00
|
|
|
|
|
|
|
// now construct a matrix with prescribed determinant
|
|
|
|
m1.setZero();
|
|
|
|
for(int i = 0; i < size; i++) m1(i,i) = ei_random<Scalar>();
|
|
|
|
RealScalar absdet = ei_abs(m1.diagonal().prod());
|
|
|
|
m3 = qr.matrixQ(); // get a unitary
|
|
|
|
m1 = m3 * m1 * m3;
|
|
|
|
qr.compute(m1);
|
|
|
|
VERIFY_IS_APPROX(absdet, qr.absDeterminant());
|
|
|
|
VERIFY_IS_APPROX(ei_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
|
|
|
ColPivHouseholderQR<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-25 01:46:14 +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())
|
|
|
|
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_colpivoting()
|
2009-08-22 13:13:21 +08:00
|
|
|
{
|
2009-12-02 02:26:29 +08:00
|
|
|
for(int i = 0; i < g_repeat; i++) {
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_1( qr<MatrixXf>() );
|
|
|
|
CALL_SUBTEST_2( qr<MatrixXd>() );
|
|
|
|
CALL_SUBTEST_3( qr<MatrixXcd>() );
|
|
|
|
CALL_SUBTEST_4(( qr_fixedsize<Matrix<float,3,5>, 4 >() ));
|
|
|
|
CALL_SUBTEST_5(( qr_fixedsize<Matrix<double,6,2>, 3 >() ));
|
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_6( 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_7(qr_verify_assert<Matrix3f>());
|
|
|
|
CALL_SUBTEST_8(qr_verify_assert<Matrix3d>());
|
|
|
|
CALL_SUBTEST_1(qr_verify_assert<MatrixXf>());
|
|
|
|
CALL_SUBTEST_2(qr_verify_assert<MatrixXd>());
|
|
|
|
CALL_SUBTEST_6(qr_verify_assert<MatrixXcf>());
|
|
|
|
CALL_SUBTEST_3(qr_verify_assert<MatrixXcd>());
|
2009-08-22 13:13:21 +08:00
|
|
|
}
|