2008-08-10 03:26:14 +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-08-10 03:26:14 +08:00
|
|
|
//
|
2009-09-22 12:16:51 +08:00
|
|
|
// Copyright (C) 2008-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
|
2008-08-10 03:26:14 +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/LU>
|
2009-11-17 06:05:12 +08:00
|
|
|
using namespace std;
|
2008-08-10 03:26:14 +08:00
|
|
|
|
|
|
|
template<typename MatrixType> void lu_non_invertible()
|
|
|
|
{
|
2009-11-17 06:05:12 +08:00
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
2008-08-10 03:26:14 +08:00
|
|
|
/* this test covers the following files:
|
|
|
|
LU.h
|
|
|
|
*/
|
2009-10-19 22:56:37 +08:00
|
|
|
int rows, cols, cols2;
|
|
|
|
if(MatrixType::RowsAtCompileTime==Dynamic)
|
|
|
|
{
|
|
|
|
rows = ei_random<int>(20,200);
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
rows = MatrixType::RowsAtCompileTime;
|
|
|
|
}
|
|
|
|
if(MatrixType::ColsAtCompileTime==Dynamic)
|
|
|
|
{
|
|
|
|
cols = ei_random<int>(20,200);
|
|
|
|
cols2 = ei_random<int>(20,200);
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
cols2 = cols = MatrixType::ColsAtCompileTime;
|
|
|
|
}
|
|
|
|
|
2009-11-09 20:51:31 +08:00
|
|
|
typedef typename ei_kernel_retval_base<FullPivLU<MatrixType> >::ReturnMatrixType KernelMatrixType;
|
|
|
|
typedef typename ei_image_retval_base<FullPivLU<MatrixType> >::ReturnMatrixType ImageMatrixType;
|
2009-10-19 22:56:37 +08:00
|
|
|
typedef Matrix<typename MatrixType::Scalar, Dynamic, Dynamic> DynamicMatrixType;
|
|
|
|
typedef Matrix<typename MatrixType::Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime>
|
|
|
|
CMatrixType;
|
|
|
|
|
2008-08-10 03:26:14 +08:00
|
|
|
int rank = ei_random<int>(1, std::min(rows, cols)-1);
|
|
|
|
|
2009-10-19 22:56:37 +08:00
|
|
|
// The image of the zero matrix should consist of a single (zero) column vector
|
2009-10-29 06:19:29 +08:00
|
|
|
VERIFY((MatrixType::Zero(rows,cols).fullPivLu().image(MatrixType::Zero(rows,cols)).cols() == 1));
|
2009-10-19 22:56:37 +08:00
|
|
|
|
|
|
|
MatrixType m1(rows, cols), m3(rows, cols2);
|
|
|
|
CMatrixType m2(cols, cols2);
|
2009-05-17 22:07:12 +08:00
|
|
|
createRandomMatrixOfRank(rank, rows, cols, m1);
|
2008-08-10 03:26:14 +08:00
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
FullPivLU<MatrixType> lu(m1);
|
2009-11-17 06:05:12 +08:00
|
|
|
// FIXME need better way to construct trapezoid matrices. extend triangularView to support rectangular.
|
|
|
|
DynamicMatrixType u(rows,cols);
|
|
|
|
for(int i = 0; i < rows; i++)
|
|
|
|
for(int j = 0; j < cols; j++)
|
|
|
|
u(i,j) = i>j ? Scalar(0) : lu.matrixLU()(i,j);
|
|
|
|
DynamicMatrixType l(rows,rows);
|
|
|
|
for(int i = 0; i < rows; i++)
|
|
|
|
for(int j = 0; j < rows; j++)
|
|
|
|
l(i,j) = (i<j || j>=cols) ? Scalar(0)
|
|
|
|
: i==j ? Scalar(1)
|
|
|
|
: lu.matrixLU()(i,j);
|
|
|
|
|
|
|
|
VERIFY_IS_APPROX(lu.permutationP() * m1 * lu.permutationQ(), l*u);
|
|
|
|
|
2009-10-19 22:56:37 +08:00
|
|
|
KernelMatrixType m1kernel = lu.kernel();
|
2009-11-03 15:18:10 +08:00
|
|
|
ImageMatrixType m1image = lu.image(m1);
|
2008-12-18 00:47:55 +08:00
|
|
|
|
2008-08-10 03:26:14 +08:00
|
|
|
VERIFY(rank == lu.rank());
|
2009-01-21 03:06:39 +08:00
|
|
|
VERIFY(cols - lu.rank() == lu.dimensionOfKernel());
|
2008-08-10 03:26:14 +08:00
|
|
|
VERIFY(!lu.isInjective());
|
|
|
|
VERIFY(!lu.isInvertible());
|
2009-08-24 12:09:01 +08:00
|
|
|
VERIFY(!lu.isSurjective());
|
2008-12-18 00:47:55 +08:00
|
|
|
VERIFY((m1 * m1kernel).isMuchSmallerThan(m1));
|
2009-10-29 06:19:29 +08:00
|
|
|
VERIFY(m1image.fullPivLu().rank() == rank);
|
2009-10-19 22:56:37 +08:00
|
|
|
DynamicMatrixType sidebyside(m1.rows(), m1.cols() + m1image.cols());
|
2008-12-18 00:47:55 +08:00
|
|
|
sidebyside << m1, m1image;
|
2009-10-29 06:19:29 +08:00
|
|
|
VERIFY(sidebyside.fullPivLu().rank() == rank);
|
2009-10-19 22:56:37 +08:00
|
|
|
m2 = CMatrixType::Random(cols,cols2);
|
2008-08-10 03:26:14 +08:00
|
|
|
m3 = m1*m2;
|
2009-10-19 22:56:37 +08:00
|
|
|
m2 = CMatrixType::Random(cols,cols2);
|
2009-09-22 13:41:09 +08:00
|
|
|
// test that the code, which does resize(), may be applied to an xpr
|
2009-09-26 23:40:29 +08:00
|
|
|
m2.block(0,0,m2.rows(),m2.cols()) = lu.solve(m3);
|
2008-08-10 03:26:14 +08:00
|
|
|
VERIFY_IS_APPROX(m3, m1*m2);
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename MatrixType> void lu_invertible()
|
|
|
|
{
|
|
|
|
/* this test covers the following files:
|
|
|
|
LU.h
|
|
|
|
*/
|
2008-08-23 01:48:36 +08:00
|
|
|
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
|
2008-08-10 03:26:14 +08:00
|
|
|
int size = ei_random<int>(10,200);
|
|
|
|
|
|
|
|
MatrixType m1(size, size), m2(size, size), m3(size, size);
|
2008-09-02 01:31:21 +08:00
|
|
|
m1 = MatrixType::Random(size,size);
|
2008-08-10 03:26:14 +08:00
|
|
|
|
2008-08-23 23:14:20 +08:00
|
|
|
if (ei_is_same_type<RealScalar,float>::ret)
|
|
|
|
{
|
|
|
|
// let's build a matrix more stable to inverse
|
2008-09-02 01:31:21 +08:00
|
|
|
MatrixType a = MatrixType::Random(size,size*2);
|
2008-08-23 23:14:20 +08:00
|
|
|
m1 += a * a.adjoint();
|
|
|
|
}
|
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
FullPivLU<MatrixType> lu(m1);
|
2008-08-10 03:26:14 +08:00
|
|
|
VERIFY(0 == lu.dimensionOfKernel());
|
2009-10-19 22:56:37 +08:00
|
|
|
VERIFY(lu.kernel().cols() == 1); // the kernel() should consist of a single (zero) column vector
|
2008-08-10 03:26:14 +08:00
|
|
|
VERIFY(size == lu.rank());
|
|
|
|
VERIFY(lu.isInjective());
|
|
|
|
VERIFY(lu.isSurjective());
|
|
|
|
VERIFY(lu.isInvertible());
|
2009-10-29 06:19:29 +08:00
|
|
|
VERIFY(lu.image(m1).fullPivLu().isInvertible());
|
2008-09-02 01:31:21 +08:00
|
|
|
m3 = MatrixType::Random(size,size);
|
2009-09-23 08:58:29 +08:00
|
|
|
m2 = lu.solve(m3);
|
2008-08-23 23:14:20 +08:00
|
|
|
VERIFY_IS_APPROX(m3, m1*m2);
|
2008-08-10 03:26:14 +08:00
|
|
|
VERIFY_IS_APPROX(m2, lu.inverse()*m3);
|
|
|
|
}
|
|
|
|
|
2009-05-22 20:27:58 +08:00
|
|
|
template<typename MatrixType> void lu_verify_assert()
|
|
|
|
{
|
|
|
|
MatrixType tmp;
|
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
FullPivLU<MatrixType> lu;
|
2009-05-22 20:27:58 +08:00
|
|
|
VERIFY_RAISES_ASSERT(lu.matrixLU())
|
|
|
|
VERIFY_RAISES_ASSERT(lu.permutationP())
|
|
|
|
VERIFY_RAISES_ASSERT(lu.permutationQ())
|
|
|
|
VERIFY_RAISES_ASSERT(lu.kernel())
|
2009-10-19 03:21:19 +08:00
|
|
|
VERIFY_RAISES_ASSERT(lu.image(tmp))
|
2009-09-23 08:58:29 +08:00
|
|
|
VERIFY_RAISES_ASSERT(lu.solve(tmp))
|
2009-05-22 20:27:58 +08:00
|
|
|
VERIFY_RAISES_ASSERT(lu.determinant())
|
|
|
|
VERIFY_RAISES_ASSERT(lu.rank())
|
|
|
|
VERIFY_RAISES_ASSERT(lu.dimensionOfKernel())
|
|
|
|
VERIFY_RAISES_ASSERT(lu.isInjective())
|
|
|
|
VERIFY_RAISES_ASSERT(lu.isSurjective())
|
|
|
|
VERIFY_RAISES_ASSERT(lu.isInvertible())
|
|
|
|
VERIFY_RAISES_ASSERT(lu.inverse())
|
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
PartialPivLU<MatrixType> plu;
|
2009-05-22 20:27:58 +08:00
|
|
|
VERIFY_RAISES_ASSERT(plu.matrixLU())
|
|
|
|
VERIFY_RAISES_ASSERT(plu.permutationP())
|
2009-10-29 06:19:29 +08:00
|
|
|
VERIFY_RAISES_ASSERT(plu.solve(tmp))
|
2009-05-22 20:27:58 +08:00
|
|
|
VERIFY_RAISES_ASSERT(plu.determinant())
|
|
|
|
VERIFY_RAISES_ASSERT(plu.inverse())
|
|
|
|
}
|
|
|
|
|
2008-08-10 03:26:14 +08:00
|
|
|
void test_lu()
|
|
|
|
{
|
|
|
|
for(int i = 0; i < g_repeat; i++) {
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_1( lu_non_invertible<Matrix3f>() );
|
|
|
|
CALL_SUBTEST_1( lu_verify_assert<Matrix3f>() );
|
2009-10-20 02:40:35 +08:00
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_2( (lu_non_invertible<Matrix<double, 4, 6> >()) );
|
|
|
|
CALL_SUBTEST_2( (lu_verify_assert<Matrix<double, 4, 6> >()) );
|
2009-10-20 02:40:35 +08:00
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_3( lu_non_invertible<MatrixXf>() );
|
|
|
|
CALL_SUBTEST_3( lu_invertible<MatrixXf>() );
|
|
|
|
CALL_SUBTEST_3( lu_verify_assert<MatrixXf>() );
|
2009-10-20 02:40:35 +08:00
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_4( lu_non_invertible<MatrixXd>() );
|
|
|
|
CALL_SUBTEST_4( lu_invertible<MatrixXd>() );
|
|
|
|
CALL_SUBTEST_4( lu_verify_assert<MatrixXd>() );
|
2009-10-20 02:40:35 +08:00
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_5( lu_non_invertible<MatrixXcf>() );
|
|
|
|
CALL_SUBTEST_5( lu_invertible<MatrixXcf>() );
|
|
|
|
CALL_SUBTEST_5( lu_verify_assert<MatrixXcf>() );
|
2009-10-20 02:40:35 +08:00
|
|
|
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_6( lu_non_invertible<MatrixXcd>() );
|
|
|
|
CALL_SUBTEST_6( lu_invertible<MatrixXcd>() );
|
|
|
|
CALL_SUBTEST_6( lu_verify_assert<MatrixXcd>() );
|
2008-08-10 03:26:14 +08:00
|
|
|
}
|
|
|
|
}
|