test/lu.cpp
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #include "main.h"
11 #include <Eigen/LU>
12 using namespace std;
13 
14 template<typename MatrixType>
16  return m.cwiseAbs().colwise().sum().maxCoeff();
17 }
18 
19 template<typename MatrixType> void lu_non_invertible()
20 {
21  typedef typename MatrixType::RealScalar RealScalar;
22  /* this test covers the following files:
23  LU.h
24  */
25  Index rows, cols, cols2;
26  if(MatrixType::RowsAtCompileTime==Dynamic)
27  {
28  rows = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE);
29  }
30  else
31  {
32  rows = MatrixType::RowsAtCompileTime;
33  }
34  if(MatrixType::ColsAtCompileTime==Dynamic)
35  {
36  cols = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE);
37  cols2 = internal::random<int>(2,EIGEN_TEST_MAX_SIZE);
38  }
39  else
40  {
41  cols2 = cols = MatrixType::ColsAtCompileTime;
42  }
43 
44  enum {
45  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
46  ColsAtCompileTime = MatrixType::ColsAtCompileTime
47  };
48  typedef typename internal::kernel_retval_base<FullPivLU<MatrixType> >::ReturnType KernelMatrixType;
49  typedef typename internal::image_retval_base<FullPivLU<MatrixType> >::ReturnType ImageMatrixType;
51  CMatrixType;
53  RMatrixType;
54 
55  Index rank = internal::random<Index>(1, (std::min)(rows, cols)-1);
56 
57  // The image of the zero matrix should consist of a single (zero) column vector
58  VERIFY((MatrixType::Zero(rows,cols).fullPivLu().image(MatrixType::Zero(rows,cols)).cols() == 1));
59 
60  // The kernel of the zero matrix is the entire space, and thus is an invertible matrix of dimensions cols.
61  KernelMatrixType kernel = MatrixType::Zero(rows,cols).fullPivLu().kernel();
62  VERIFY((kernel.fullPivLu().isInvertible()));
63 
64  MatrixType m1(rows, cols), m3(rows, cols2);
65  CMatrixType m2(cols, cols2);
66  createRandomPIMatrixOfRank(rank, rows, cols, m1);
67 
69 
70  // The special value 0.01 below works well in tests. Keep in mind that we're only computing the rank
71  // of singular values are either 0 or 1.
72  // So it's not clear at all that the epsilon should play any role there.
73  lu.setThreshold(RealScalar(0.01));
74  lu.compute(m1);
75 
76  MatrixType u(rows,cols);
77  u = lu.matrixLU().template triangularView<Upper>();
78  RMatrixType l = RMatrixType::Identity(rows,rows);
79  l.block(0,0,rows,(std::min)(rows,cols)).template triangularView<StrictlyLower>()
80  = lu.matrixLU().block(0,0,rows,(std::min)(rows,cols));
81 
82  VERIFY_IS_APPROX(lu.permutationP() * m1 * lu.permutationQ(), l*u);
83 
84  KernelMatrixType m1kernel = lu.kernel();
85  ImageMatrixType m1image = lu.image(m1);
86 
88  VERIFY(rank == lu.rank());
89  VERIFY(cols - lu.rank() == lu.dimensionOfKernel());
90  VERIFY(!lu.isInjective());
91  VERIFY(!lu.isInvertible());
92  VERIFY(!lu.isSurjective());
93  VERIFY_IS_MUCH_SMALLER_THAN((m1 * m1kernel), m1);
94  VERIFY(m1image.fullPivLu().rank() == rank);
95  VERIFY_IS_APPROX(m1 * m1.adjoint() * m1image, m1image);
96 
97  m2 = CMatrixType::Random(cols,cols2);
98  m3 = m1*m2;
99  m2 = CMatrixType::Random(cols,cols2);
100  // test that the code, which does resize(), may be applied to an xpr
101  m2.block(0,0,m2.rows(),m2.cols()) = lu.solve(m3);
102  VERIFY_IS_APPROX(m3, m1*m2);
103 
104  // test solve with transposed
105  m3 = MatrixType::Random(rows,cols2);
106  m2 = m1.transpose()*m3;
107  m3 = MatrixType::Random(rows,cols2);
108  lu.template _solve_impl_transposed<false>(m2, m3);
109  VERIFY_IS_APPROX(m2, m1.transpose()*m3);
110  m3 = MatrixType::Random(rows,cols2);
111  m3 = lu.transpose().solve(m2);
112  VERIFY_IS_APPROX(m2, m1.transpose()*m3);
113 
114  // test solve with conjugate transposed
115  m3 = MatrixType::Random(rows,cols2);
116  m2 = m1.adjoint()*m3;
117  m3 = MatrixType::Random(rows,cols2);
118  lu.template _solve_impl_transposed<true>(m2, m3);
119  VERIFY_IS_APPROX(m2, m1.adjoint()*m3);
120  m3 = MatrixType::Random(rows,cols2);
121  m3 = lu.adjoint().solve(m2);
122  VERIFY_IS_APPROX(m2, m1.adjoint()*m3);
123 }
124 
125 template<typename MatrixType> void lu_invertible()
126 {
127  /* this test covers the following files:
128  LU.h
129  */
131  Index size = MatrixType::RowsAtCompileTime;
132  if( size==Dynamic)
133  size = internal::random<Index>(1,EIGEN_TEST_MAX_SIZE);
134 
135  MatrixType m1(size, size), m2(size, size), m3(size, size);
137  lu.setThreshold(RealScalar(0.01));
138  do {
139  m1 = MatrixType::Random(size,size);
140  lu.compute(m1);
141  } while(!lu.isInvertible());
142 
144  VERIFY(0 == lu.dimensionOfKernel());
145  VERIFY(lu.kernel().cols() == 1); // the kernel() should consist of a single (zero) column vector
146  VERIFY(size == lu.rank());
147  VERIFY(lu.isInjective());
148  VERIFY(lu.isSurjective());
149  VERIFY(lu.isInvertible());
150  VERIFY(lu.image(m1).fullPivLu().isInvertible());
151  m3 = MatrixType::Random(size,size);
152  m2 = lu.solve(m3);
153  VERIFY_IS_APPROX(m3, m1*m2);
154  MatrixType m1_inverse = lu.inverse();
155  VERIFY_IS_APPROX(m2, m1_inverse*m3);
156 
157  RealScalar rcond = (RealScalar(1) / matrix_l1_norm(m1)) / matrix_l1_norm(m1_inverse);
158  const RealScalar rcond_est = lu.rcond();
159  // Verify that the estimated condition number is within a factor of 10 of the
160  // truth.
161  VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10);
162 
163  // test solve with transposed
164  lu.template _solve_impl_transposed<false>(m3, m2);
165  VERIFY_IS_APPROX(m3, m1.transpose()*m2);
166  m3 = MatrixType::Random(size,size);
167  m3 = lu.transpose().solve(m2);
168  VERIFY_IS_APPROX(m2, m1.transpose()*m3);
169 
170  // test solve with conjugate transposed
171  lu.template _solve_impl_transposed<true>(m3, m2);
172  VERIFY_IS_APPROX(m3, m1.adjoint()*m2);
173  m3 = MatrixType::Random(size,size);
174  m3 = lu.adjoint().solve(m2);
175  VERIFY_IS_APPROX(m2, m1.adjoint()*m3);
176 
177  // Regression test for Bug 302
178  MatrixType m4 = MatrixType::Random(size,size);
179  VERIFY_IS_APPROX(lu.solve(m3*m4), lu.solve(m3)*m4);
180 }
181 
182 template<typename MatrixType> void lu_partial_piv()
183 {
184  /* this test covers the following files:
185  PartialPivLU.h
186  */
188  Index size = internal::random<Index>(1,4);
189 
190  MatrixType m1(size, size), m2(size, size), m3(size, size);
191  m1.setRandom();
193 
195 
196  m3 = MatrixType::Random(size,size);
197  m2 = plu.solve(m3);
198  VERIFY_IS_APPROX(m3, m1*m2);
199  MatrixType m1_inverse = plu.inverse();
200  VERIFY_IS_APPROX(m2, m1_inverse*m3);
201 
202  RealScalar rcond = (RealScalar(1) / matrix_l1_norm(m1)) / matrix_l1_norm(m1_inverse);
203  const RealScalar rcond_est = plu.rcond();
204  // Verify that the estimate is within a factor of 10 of the truth.
205  VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10);
206 
207  // test solve with transposed
208  plu.template _solve_impl_transposed<false>(m3, m2);
209  VERIFY_IS_APPROX(m3, m1.transpose()*m2);
210  m3 = MatrixType::Random(size,size);
211  m3 = plu.transpose().solve(m2);
212  VERIFY_IS_APPROX(m2, m1.transpose()*m3);
213 
214  // test solve with conjugate transposed
215  plu.template _solve_impl_transposed<true>(m3, m2);
216  VERIFY_IS_APPROX(m3, m1.adjoint()*m2);
217  m3 = MatrixType::Random(size,size);
218  m3 = plu.adjoint().solve(m2);
219  VERIFY_IS_APPROX(m2, m1.adjoint()*m3);
220 }
221 
222 template<typename MatrixType> void lu_verify_assert()
223 {
224  MatrixType tmp;
225 
231  VERIFY_RAISES_ASSERT(lu.image(tmp))
232  VERIFY_RAISES_ASSERT(lu.solve(tmp))
240 
242  VERIFY_RAISES_ASSERT(plu.matrixLU())
243  VERIFY_RAISES_ASSERT(plu.permutationP())
244  VERIFY_RAISES_ASSERT(plu.solve(tmp))
245  VERIFY_RAISES_ASSERT(plu.determinant())
246  VERIFY_RAISES_ASSERT(plu.inverse())
247 }
248 
249 void test_lu()
250 {
251  for(int i = 0; i < g_repeat; i++) {
252  CALL_SUBTEST_1( lu_non_invertible<Matrix3f>() );
253  CALL_SUBTEST_1( lu_invertible<Matrix3f>() );
254  CALL_SUBTEST_1( lu_verify_assert<Matrix3f>() );
255 
256  CALL_SUBTEST_2( (lu_non_invertible<Matrix<double, 4, 6> >()) );
257  CALL_SUBTEST_2( (lu_verify_assert<Matrix<double, 4, 6> >()) );
258 
259  CALL_SUBTEST_3( lu_non_invertible<MatrixXf>() );
260  CALL_SUBTEST_3( lu_invertible<MatrixXf>() );
261  CALL_SUBTEST_3( lu_verify_assert<MatrixXf>() );
262 
263  CALL_SUBTEST_4( lu_non_invertible<MatrixXd>() );
264  CALL_SUBTEST_4( lu_invertible<MatrixXd>() );
265  CALL_SUBTEST_4( lu_partial_piv<MatrixXd>() );
266  CALL_SUBTEST_4( lu_verify_assert<MatrixXd>() );
267 
268  CALL_SUBTEST_5( lu_non_invertible<MatrixXcf>() );
269  CALL_SUBTEST_5( lu_invertible<MatrixXcf>() );
270  CALL_SUBTEST_5( lu_verify_assert<MatrixXcf>() );
271 
272  CALL_SUBTEST_6( lu_non_invertible<MatrixXcd>() );
273  CALL_SUBTEST_6( lu_invertible<MatrixXcd>() );
274  CALL_SUBTEST_6( lu_partial_piv<MatrixXcd>() );
275  CALL_SUBTEST_6( lu_verify_assert<MatrixXcd>() );
276 
277  CALL_SUBTEST_7(( lu_non_invertible<Matrix<float,Dynamic,16> >() ));
278 
279  // Test problem size constructors
280  CALL_SUBTEST_9( PartialPivLU<MatrixXf>(10) );
281  CALL_SUBTEST_9( FullPivLU<MatrixXf>(10, 20); );
282  }
283 }
Matrix3f m
void lu_invertible()
Definition: test/lu.cpp:125
#define VERIFY_RAISES_ASSERT(a)
Definition: main.h:285
bool isInjective() const
Definition: FullPivLU.h:362
RealScalar rcond() const
Definition: PartialPivLU.h:184
const Inverse< PartialPivLU > inverse() const
Definition: PartialPivLU.h:197
const internal::image_retval< FullPivLU > image(const MatrixType &originalMatrix) const
Definition: FullPivLU.h:215
FullPivLU & setThreshold(const RealScalar &threshold)
Definition: FullPivLU.h:292
internal::traits< MatrixType >::Scalar determinant() const
Definition: FullPivLU.h:583
void lu_verify_assert()
Definition: test/lu.cpp:222
#define min(a, b)
Definition: datatypes.h:19
RealScalar rcond() const
Definition: FullPivLU.h:252
MatrixType m2(n_dims)
EIGEN_DEVICE_FUNC const PermutationPType & permutationP() const
Definition: FullPivLU.h:159
FullPivLU & compute(const EigenBase< InputType > &matrix)
Definition: FullPivLU.h:119
LU decomposition of a matrix with partial pivoting, and related features.
Definition: Half.h:150
MatrixXf MatrixType
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:150
void lu_non_invertible()
Definition: test/lu.cpp:19
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
bool isSurjective() const
Definition: FullPivLU.h:375
#define VERIFY_IS_APPROX(a, b)
static const Line3 l(Rot3(), 1, 1)
Matrix3d m1
Definition: IOFormat.cpp:2
static int g_repeat
Definition: main.h:144
const PermutationQType & permutationQ() const
Definition: FullPivLU.h:169
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
cout<< "Here is the matrix m:"<< endl<< m<< endl;Eigen::FullPivLU< Matrix5x3 > lu(m)
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:34
MatrixType reconstructedMatrix() const
Definition: FullPivLU.h:594
#define VERIFY_IS_MUCH_SMALLER_THAN(a, b)
Definition: main.h:335
ConstTransposeReturnType transpose() const
Definition: SolverBase.h:90
void test_lu()
Definition: test/lu.cpp:249
#define VERIFY(a)
Definition: main.h:325
#define EIGEN_TEST_MAX_SIZE
LU decomposition of a matrix with complete pivoting, and related features.
bool isInvertible() const
Definition: FullPivLU.h:387
MatrixType::RealScalar matrix_l1_norm(const MatrixType &m)
Definition: test/lu.cpp:15
Index rank() const
Definition: FullPivLU.h:332
void createRandomPIMatrixOfRank(Index desired_rank, Index rows, Index cols, MatrixType &m)
Definition: main.h:603
void lu_partial_piv()
Definition: test/lu.cpp:182
const MatrixType & matrixLU() const
Definition: FullPivLU.h:131
MatrixType reconstructedMatrix() const
Definition: PartialPivLU.h:549
const int Dynamic
Definition: Constants.h:21
The matrix class, also used for vectors and row-vectors.
const Inverse< FullPivLU > inverse() const
Definition: FullPivLU.h:400
const Solve< PartialPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: PartialPivLU.h:175
Index dimensionOfKernel() const
Definition: FullPivLU.h:349
const Solve< FullPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: FullPivLU.h:243
const internal::kernel_retval< FullPivLU > kernel() const
Definition: FullPivLU.h:189


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:42:35