2008-08-09 19:26:14 +00:00
|
|
|
// This file is part of Eigen, a lightweight C++ template library
|
2009-05-22 20:25:33 +02:00
|
|
|
// for linear algebra.
|
2008-08-09 19:26:14 +00:00
|
|
|
//
|
2009-09-22 00:16:51 -04:00
|
|
|
// Copyright (C) 2008-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
|
2008-08-09 19:26:14 +00: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-16 17:05:12 -05:00
|
|
|
using namespace std;
|
2008-08-09 19:26:14 +00:00
|
|
|
|
|
|
|
|
template<typename MatrixType> void lu_non_invertible()
|
|
|
|
|
{
|
2009-11-16 17:05:12 -05:00
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
2010-02-23 09:04:59 -05:00
|
|
|
typedef typename MatrixType::RealScalar RealScalar;
|
2008-08-09 19:26:14 +00:00
|
|
|
/* this test covers the following files:
|
|
|
|
|
LU.h
|
|
|
|
|
*/
|
2009-10-19 10:56:37 -04:00
|
|
|
int rows, cols, cols2;
|
|
|
|
|
if(MatrixType::RowsAtCompileTime==Dynamic)
|
|
|
|
|
{
|
2010-01-27 11:42:04 -05:00
|
|
|
rows = ei_random<int>(2,200);
|
2009-10-19 10:56:37 -04:00
|
|
|
}
|
|
|
|
|
else
|
|
|
|
|
{
|
|
|
|
|
rows = MatrixType::RowsAtCompileTime;
|
|
|
|
|
}
|
|
|
|
|
if(MatrixType::ColsAtCompileTime==Dynamic)
|
|
|
|
|
{
|
2010-01-27 11:42:04 -05:00
|
|
|
cols = ei_random<int>(2,200);
|
|
|
|
|
cols2 = ei_random<int>(2,200);
|
2009-10-19 10:56:37 -04:00
|
|
|
}
|
|
|
|
|
else
|
|
|
|
|
{
|
|
|
|
|
cols2 = cols = MatrixType::ColsAtCompileTime;
|
|
|
|
|
}
|
|
|
|
|
|
2010-02-23 15:40:24 -05:00
|
|
|
enum {
|
|
|
|
|
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
|
|
|
|
ColsAtCompileTime = MatrixType::ColsAtCompileTime
|
|
|
|
|
};
|
2010-02-20 15:53:57 +01:00
|
|
|
typedef typename ei_kernel_retval_base<FullPivLU<MatrixType> >::ReturnType KernelMatrixType;
|
|
|
|
|
typedef typename ei_image_retval_base<FullPivLU<MatrixType> >::ReturnType ImageMatrixType;
|
2010-02-23 15:40:24 -05:00
|
|
|
typedef Matrix<typename MatrixType::Scalar, ColsAtCompileTime, ColsAtCompileTime>
|
2009-10-19 10:56:37 -04:00
|
|
|
CMatrixType;
|
2010-02-23 15:40:24 -05:00
|
|
|
typedef Matrix<typename MatrixType::Scalar, RowsAtCompileTime, RowsAtCompileTime>
|
|
|
|
|
RMatrixType;
|
2009-10-19 10:56:37 -04:00
|
|
|
|
2008-08-09 19:26:14 +00:00
|
|
|
int rank = ei_random<int>(1, std::min(rows, cols)-1);
|
|
|
|
|
|
2009-10-19 10:56:37 -04:00
|
|
|
// The image of the zero matrix should consist of a single (zero) column vector
|
2009-10-28 18:19:29 -04:00
|
|
|
VERIFY((MatrixType::Zero(rows,cols).fullPivLu().image(MatrixType::Zero(rows,cols)).cols() == 1));
|
2010-03-03 18:47:58 +01:00
|
|
|
|
2009-10-19 10:56:37 -04:00
|
|
|
MatrixType m1(rows, cols), m3(rows, cols2);
|
|
|
|
|
CMatrixType m2(cols, cols2);
|
2010-02-23 15:40:24 -05:00
|
|
|
createRandomPIMatrixOfRank(rank, rows, cols, m1);
|
2010-02-23 09:04:59 -05:00
|
|
|
|
|
|
|
|
FullPivLU<MatrixType> lu;
|
|
|
|
|
|
2010-02-23 15:40:24 -05:00
|
|
|
// The special value 0.01 below works well in tests. Keep in mind that we're only computing the rank
|
|
|
|
|
// of singular values are either 0 or 1.
|
|
|
|
|
// So it's not clear at all that the epsilon should play any role there.
|
2010-02-23 09:04:59 -05:00
|
|
|
lu.setThreshold(RealScalar(0.01));
|
|
|
|
|
lu.compute(m1);
|
2008-08-09 19:26:14 +00:00
|
|
|
|
2010-02-23 15:40:24 -05:00
|
|
|
MatrixType u(rows,cols);
|
|
|
|
|
u = lu.matrixLU().template triangularView<Upper>();
|
|
|
|
|
RMatrixType l = RMatrixType::Identity(rows,rows);
|
|
|
|
|
l.block(0,0,rows,std::min(rows,cols)).template triangularView<StrictlyLower>()
|
|
|
|
|
= lu.matrixLU().block(0,0,rows,std::min(rows,cols));
|
2010-03-03 18:47:58 +01:00
|
|
|
|
2009-11-16 17:05:12 -05:00
|
|
|
VERIFY_IS_APPROX(lu.permutationP() * m1 * lu.permutationQ(), l*u);
|
2010-03-03 18:47:58 +01:00
|
|
|
|
2009-10-19 10:56:37 -04:00
|
|
|
KernelMatrixType m1kernel = lu.kernel();
|
2009-11-03 02:18:10 -05:00
|
|
|
ImageMatrixType m1image = lu.image(m1);
|
2008-12-17 16:47:55 +00:00
|
|
|
|
2010-02-24 19:16:10 +01:00
|
|
|
VERIFY_IS_APPROX(m1, lu.reconstructedMatrix());
|
2008-08-09 19:26:14 +00:00
|
|
|
VERIFY(rank == lu.rank());
|
2009-01-20 19:06:39 +00:00
|
|
|
VERIFY(cols - lu.rank() == lu.dimensionOfKernel());
|
2008-08-09 19:26:14 +00:00
|
|
|
VERIFY(!lu.isInjective());
|
|
|
|
|
VERIFY(!lu.isInvertible());
|
2009-08-24 00:09:01 -04:00
|
|
|
VERIFY(!lu.isSurjective());
|
2008-12-17 16:47:55 +00:00
|
|
|
VERIFY((m1 * m1kernel).isMuchSmallerThan(m1));
|
2009-10-28 18:19:29 -04:00
|
|
|
VERIFY(m1image.fullPivLu().rank() == rank);
|
2010-02-23 15:40:24 -05:00
|
|
|
VERIFY_IS_APPROX(m1 * m1.adjoint() * m1image, m1image);
|
2010-02-23 09:04:59 -05:00
|
|
|
|
2009-10-19 10:56:37 -04:00
|
|
|
m2 = CMatrixType::Random(cols,cols2);
|
2008-08-09 19:26:14 +00:00
|
|
|
m3 = m1*m2;
|
2009-10-19 10:56:37 -04:00
|
|
|
m2 = CMatrixType::Random(cols,cols2);
|
2009-09-22 01:41:09 -04:00
|
|
|
// test that the code, which does resize(), may be applied to an xpr
|
2009-09-26 11:40:29 -04:00
|
|
|
m2.block(0,0,m2.rows(),m2.cols()) = lu.solve(m3);
|
2008-08-09 19:26:14 +00:00
|
|
|
VERIFY_IS_APPROX(m3, m1*m2);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
template<typename MatrixType> void lu_invertible()
|
|
|
|
|
{
|
|
|
|
|
/* this test covers the following files:
|
|
|
|
|
LU.h
|
|
|
|
|
*/
|
2010-02-23 15:40:24 -05:00
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
2008-08-22 17:48:36 +00:00
|
|
|
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
|
2010-01-27 11:42:04 -05:00
|
|
|
int size = ei_random<int>(1,200);
|
2008-08-09 19:26:14 +00:00
|
|
|
|
|
|
|
|
MatrixType m1(size, size), m2(size, size), m3(size, size);
|
2010-02-23 15:40:24 -05:00
|
|
|
FullPivLU<MatrixType> lu;
|
2010-06-08 16:02:22 +02:00
|
|
|
lu.setThreshold(RealScalar(0.01));
|
2010-02-23 15:40:24 -05:00
|
|
|
do {
|
|
|
|
|
m1 = MatrixType::Random(size,size);
|
|
|
|
|
lu.compute(m1);
|
|
|
|
|
} while(!lu.isInvertible());
|
2008-08-23 15:14:20 +00:00
|
|
|
|
2010-02-24 19:16:10 +01:00
|
|
|
VERIFY_IS_APPROX(m1, lu.reconstructedMatrix());
|
2008-08-09 19:26:14 +00:00
|
|
|
VERIFY(0 == lu.dimensionOfKernel());
|
2009-10-19 10:56:37 -04:00
|
|
|
VERIFY(lu.kernel().cols() == 1); // the kernel() should consist of a single (zero) column vector
|
2008-08-09 19:26:14 +00:00
|
|
|
VERIFY(size == lu.rank());
|
|
|
|
|
VERIFY(lu.isInjective());
|
|
|
|
|
VERIFY(lu.isSurjective());
|
|
|
|
|
VERIFY(lu.isInvertible());
|
2009-10-28 18:19:29 -04:00
|
|
|
VERIFY(lu.image(m1).fullPivLu().isInvertible());
|
2008-09-01 17:31:21 +00:00
|
|
|
m3 = MatrixType::Random(size,size);
|
2009-09-22 20:58:29 -04:00
|
|
|
m2 = lu.solve(m3);
|
2008-08-23 15:14:20 +00:00
|
|
|
VERIFY_IS_APPROX(m3, m1*m2);
|
2008-08-09 19:26:14 +00:00
|
|
|
VERIFY_IS_APPROX(m2, lu.inverse()*m3);
|
|
|
|
|
}
|
|
|
|
|
|
2010-02-24 19:16:10 +01:00
|
|
|
template<typename MatrixType> void lu_partial_piv()
|
|
|
|
|
{
|
|
|
|
|
/* this test covers the following files:
|
|
|
|
|
PartialPivLU.h
|
|
|
|
|
*/
|
|
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
|
|
|
|
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
|
|
|
|
|
int rows = ei_random<int>(1,4);
|
|
|
|
|
int cols = rows;
|
|
|
|
|
|
|
|
|
|
MatrixType m1(cols, rows);
|
|
|
|
|
m1.setRandom();
|
|
|
|
|
PartialPivLU<MatrixType> plu(m1);
|
|
|
|
|
|
|
|
|
|
VERIFY_IS_APPROX(m1, plu.reconstructedMatrix());
|
|
|
|
|
}
|
|
|
|
|
|
2009-05-22 14:27:58 +02:00
|
|
|
template<typename MatrixType> void lu_verify_assert()
|
|
|
|
|
{
|
|
|
|
|
MatrixType tmp;
|
|
|
|
|
|
2009-10-28 18:19:29 -04:00
|
|
|
FullPivLU<MatrixType> lu;
|
2009-05-22 14:27:58 +02:00
|
|
|
VERIFY_RAISES_ASSERT(lu.matrixLU())
|
|
|
|
|
VERIFY_RAISES_ASSERT(lu.permutationP())
|
|
|
|
|
VERIFY_RAISES_ASSERT(lu.permutationQ())
|
|
|
|
|
VERIFY_RAISES_ASSERT(lu.kernel())
|
2009-10-18 15:21:19 -04:00
|
|
|
VERIFY_RAISES_ASSERT(lu.image(tmp))
|
2009-09-22 20:58:29 -04:00
|
|
|
VERIFY_RAISES_ASSERT(lu.solve(tmp))
|
2009-05-22 14:27:58 +02: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-28 18:19:29 -04:00
|
|
|
PartialPivLU<MatrixType> plu;
|
2009-05-22 14:27:58 +02:00
|
|
|
VERIFY_RAISES_ASSERT(plu.matrixLU())
|
|
|
|
|
VERIFY_RAISES_ASSERT(plu.permutationP())
|
2009-10-28 18:19:29 -04:00
|
|
|
VERIFY_RAISES_ASSERT(plu.solve(tmp))
|
2009-05-22 14:27:58 +02:00
|
|
|
VERIFY_RAISES_ASSERT(plu.determinant())
|
|
|
|
|
VERIFY_RAISES_ASSERT(plu.inverse())
|
|
|
|
|
}
|
|
|
|
|
|
2008-08-09 19:26:14 +00:00
|
|
|
void test_lu()
|
|
|
|
|
{
|
|
|
|
|
for(int i = 0; i < g_repeat; i++) {
|
2009-10-28 18:19:29 -04:00
|
|
|
CALL_SUBTEST_1( lu_non_invertible<Matrix3f>() );
|
|
|
|
|
CALL_SUBTEST_1( lu_verify_assert<Matrix3f>() );
|
2009-10-19 14:40:35 -04:00
|
|
|
|
2009-10-28 18:19:29 -04:00
|
|
|
CALL_SUBTEST_2( (lu_non_invertible<Matrix<double, 4, 6> >()) );
|
|
|
|
|
CALL_SUBTEST_2( (lu_verify_assert<Matrix<double, 4, 6> >()) );
|
2010-03-03 18:47:58 +01:00
|
|
|
|
2009-10-28 18:19:29 -04:00
|
|
|
CALL_SUBTEST_3( lu_non_invertible<MatrixXf>() );
|
|
|
|
|
CALL_SUBTEST_3( lu_invertible<MatrixXf>() );
|
|
|
|
|
CALL_SUBTEST_3( lu_verify_assert<MatrixXf>() );
|
2010-03-03 18:47:58 +01:00
|
|
|
|
2009-10-28 18:19:29 -04:00
|
|
|
CALL_SUBTEST_4( lu_non_invertible<MatrixXd>() );
|
|
|
|
|
CALL_SUBTEST_4( lu_invertible<MatrixXd>() );
|
2010-02-24 19:16:10 +01:00
|
|
|
CALL_SUBTEST_4( lu_partial_piv<MatrixXd>() );
|
2009-10-28 18:19:29 -04:00
|
|
|
CALL_SUBTEST_4( lu_verify_assert<MatrixXd>() );
|
2010-03-03 18:47:58 +01:00
|
|
|
|
2009-10-28 18:19:29 -04:00
|
|
|
CALL_SUBTEST_5( lu_non_invertible<MatrixXcf>() );
|
|
|
|
|
CALL_SUBTEST_5( lu_invertible<MatrixXcf>() );
|
|
|
|
|
CALL_SUBTEST_5( lu_verify_assert<MatrixXcf>() );
|
2010-03-03 18:47:58 +01:00
|
|
|
|
2009-10-28 18:19:29 -04:00
|
|
|
CALL_SUBTEST_6( lu_non_invertible<MatrixXcd>() );
|
|
|
|
|
CALL_SUBTEST_6( lu_invertible<MatrixXcd>() );
|
2010-02-24 19:16:10 +01:00
|
|
|
CALL_SUBTEST_6( lu_partial_piv<MatrixXcd>() );
|
2009-10-28 18:19:29 -04:00
|
|
|
CALL_SUBTEST_6( lu_verify_assert<MatrixXcd>() );
|
2010-01-27 11:42:04 -05:00
|
|
|
|
|
|
|
|
CALL_SUBTEST_7(( lu_non_invertible<Matrix<float,Dynamic,16> >() ));
|
2010-04-21 17:15:57 +02:00
|
|
|
|
|
|
|
|
// Test problem size constructors
|
|
|
|
|
CALL_SUBTEST_9( PartialPivLU<MatrixXf>(10) );
|
|
|
|
|
CALL_SUBTEST_9( FullPivLU<MatrixXf>(10, 20); );
|
2008-08-09 19:26:14 +00:00
|
|
|
}
|
|
|
|
|
}
|