2008-05-27 17:16:27 +08:00
|
|
|
// This file is part of Eigen, a lightweight C++ template library
|
2009-05-23 02:25:33 +08:00
|
|
|
// for linear algebra.
|
2008-05-27 17:16:27 +08:00
|
|
|
//
|
2010-06-25 05:21:58 +08:00
|
|
|
// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
|
2008-05-27 17:16:27 +08:00
|
|
|
//
|
|
|
|
// 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(const MatrixType& m)
|
|
|
|
{
|
2010-06-20 23:37:56 +08:00
|
|
|
typedef typename MatrixType::Index Index;
|
|
|
|
|
|
|
|
Index rows = m.rows();
|
|
|
|
Index cols = m.cols();
|
2008-05-27 17:16:27 +08:00
|
|
|
|
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
2009-10-13 10:33:51 +08:00
|
|
|
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
|
2008-05-27 17:16:27 +08:00
|
|
|
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
|
|
|
|
|
2008-07-21 08:34:46 +08:00
|
|
|
MatrixType a = MatrixType::Random(rows,cols);
|
2009-07-06 23:12:10 +08:00
|
|
|
HouseholderQR<MatrixType> qrOfA(a);
|
2010-06-18 00:30:47 +08:00
|
|
|
|
2009-12-03 00:11:09 +08:00
|
|
|
MatrixQType q = qrOfA.householderQ();
|
2009-10-13 10:33:51 +08:00
|
|
|
VERIFY_IS_UNITARY(q);
|
2010-06-18 00:30:47 +08:00
|
|
|
|
2010-01-11 21:48:39 +08:00
|
|
|
MatrixType r = qrOfA.matrixQR().template triangularView<Upper>();
|
2009-12-03 00:11:09 +08:00
|
|
|
VERIFY_IS_APPROX(a, qrOfA.householderQ() * r);
|
2009-09-28 22:49:55 +08:00
|
|
|
}
|
2008-06-08 23:03:23 +08:00
|
|
|
|
2009-09-28 22:49:55 +08:00
|
|
|
template<typename MatrixType, int Cols2> void qr_fixedsize()
|
|
|
|
{
|
|
|
|
enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime };
|
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
|
|
|
Matrix<Scalar,Rows,Cols> m1 = Matrix<Scalar,Rows,Cols>::Random();
|
|
|
|
HouseholderQR<Matrix<Scalar,Rows,Cols> > qr(m1);
|
2008-06-08 23:03:23 +08:00
|
|
|
|
2009-09-28 22:49:55 +08:00
|
|
|
Matrix<Scalar,Rows,Cols> r = qr.matrixQR();
|
|
|
|
// 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);
|
2008-06-08 23:03:23 +08:00
|
|
|
|
2009-12-03 00:11:09 +08:00
|
|
|
VERIFY_IS_APPROX(m1, qr.householderQ() * r);
|
2009-09-28 22:49:55 +08:00
|
|
|
|
|
|
|
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 22:49:55 +08:00
|
|
|
VERIFY_IS_APPROX(m3, m1*m2);
|
2008-05-27 17:16:27 +08:00
|
|
|
}
|
|
|
|
|
2009-01-28 17:45:53 +08:00
|
|
|
template<typename MatrixType> void qr_invertible()
|
|
|
|
{
|
|
|
|
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
|
2009-08-25 01:46:14 +08:00
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
|
|
|
|
2010-10-25 22:15:22 +08:00
|
|
|
int size = internal::random<int>(10,50);
|
2009-01-28 17:45:53 +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-01-28 17:45:53 +08:00
|
|
|
{
|
|
|
|
// let's build a matrix more stable to inverse
|
|
|
|
MatrixType a = MatrixType::Random(size,size*2);
|
|
|
|
m1 += a * a.adjoint();
|
|
|
|
}
|
|
|
|
|
2009-07-06 23:12:10 +08:00
|
|
|
HouseholderQR<MatrixType> qr(m1);
|
2009-01-28 17:45:53 +08:00
|
|
|
m3 = MatrixType::Random(size,size);
|
2009-11-08 23:21:26 +08:00
|
|
|
m2 = qr.solve(m3);
|
2009-01-28 17:45:53 +08:00
|
|
|
VERIFY_IS_APPROX(m3, m1*m2);
|
2009-09-16 20:35:42 +08:00
|
|
|
|
2009-08-25 01:46:14 +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>();
|
|
|
|
RealScalar absdet = internal::abs(m1.diagonal().prod());
|
2009-12-03 00:11:09 +08:00
|
|
|
m3 = qr.householderQ(); // get a unitary
|
2009-08-25 01:46:14 +08:00
|
|
|
m1 = m3 * m1 * m3;
|
|
|
|
qr.compute(m1);
|
|
|
|
VERIFY_IS_APPROX(absdet, qr.absDeterminant());
|
2010-10-25 22:15:22 +08:00
|
|
|
VERIFY_IS_APPROX(internal::log(absdet), qr.logAbsDeterminant());
|
2009-01-28 17:45:53 +08:00
|
|
|
}
|
|
|
|
|
2009-05-22 20:27:58 +08:00
|
|
|
template<typename MatrixType> void qr_verify_assert()
|
|
|
|
{
|
|
|
|
MatrixType tmp;
|
|
|
|
|
2009-07-06 23:12:10 +08:00
|
|
|
HouseholderQR<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-12-03 00:11:09 +08:00
|
|
|
VERIFY_RAISES_ASSERT(qr.householderQ())
|
2009-08-25 01:46:14 +08:00
|
|
|
VERIFY_RAISES_ASSERT(qr.absDeterminant())
|
|
|
|
VERIFY_RAISES_ASSERT(qr.logAbsDeterminant())
|
2009-05-22 20:27:58 +08:00
|
|
|
}
|
|
|
|
|
2008-05-27 17:16:27 +08:00
|
|
|
void test_qr()
|
|
|
|
{
|
2010-01-11 21:48:39 +08:00
|
|
|
for(int i = 0; i < g_repeat; i++) {
|
2010-10-25 22:15:22 +08:00
|
|
|
CALL_SUBTEST_1( qr(MatrixXf(internal::random<int>(1,200),internal::random<int>(1,200))) );
|
|
|
|
CALL_SUBTEST_2( qr(MatrixXcd(internal::random<int>(1,200),internal::random<int>(1,200))) );
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_3(( qr_fixedsize<Matrix<float,3,4>, 2 >() ));
|
|
|
|
CALL_SUBTEST_4(( qr_fixedsize<Matrix<double,6,2>, 4 >() ));
|
|
|
|
CALL_SUBTEST_5(( qr_fixedsize<Matrix<double,2,5>, 7 >() ));
|
2010-03-01 21:46:41 +08:00
|
|
|
CALL_SUBTEST_11( qr(Matrix<float,1,1>()) );
|
2008-05-27 17:16:27 +08:00
|
|
|
}
|
2009-02-04 23:37:00 +08:00
|
|
|
|
2009-01-28 17:45:53 +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_6( qr_invertible<MatrixXd>() );
|
|
|
|
CALL_SUBTEST_7( qr_invertible<MatrixXcf>() );
|
|
|
|
CALL_SUBTEST_8( qr_invertible<MatrixXcd>() );
|
2009-01-20 18:38:56 +08:00
|
|
|
}
|
2009-05-22 20:27:58 +08:00
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_9(qr_verify_assert<Matrix3f>());
|
|
|
|
CALL_SUBTEST_10(qr_verify_assert<Matrix3d>());
|
|
|
|
CALL_SUBTEST_1(qr_verify_assert<MatrixXf>());
|
|
|
|
CALL_SUBTEST_6(qr_verify_assert<MatrixXd>());
|
|
|
|
CALL_SUBTEST_7(qr_verify_assert<MatrixXcf>());
|
|
|
|
CALL_SUBTEST_8(qr_verify_assert<MatrixXcd>());
|
2010-04-21 23:15:57 +08:00
|
|
|
|
|
|
|
// Test problem size constructors
|
|
|
|
CALL_SUBTEST_12(HouseholderQR<MatrixXf>(10, 20));
|
2008-05-27 17:16:27 +08:00
|
|
|
}
|