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 #include "solverbase.h"
13 using namespace std;
14 
15 template<typename MatrixType>
17  return m.cwiseAbs().colwise().sum().maxCoeff();
18 }
19 
20 template<typename MatrixType> void lu_non_invertible()
21 {
22  STATIC_CHECK(( internal::is_same<typename FullPivLU<MatrixType>::StorageIndex,int>::value ));
23 
24  typedef typename MatrixType::RealScalar RealScalar;
25  /* this test covers the following files:
26  LU.h
27  */
28  Index rows, cols, cols2;
29  if(MatrixType::RowsAtCompileTime==Dynamic)
30  {
31  rows = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE);
32  }
33  else
34  {
35  rows = MatrixType::RowsAtCompileTime;
36  }
37  if(MatrixType::ColsAtCompileTime==Dynamic)
38  {
39  cols = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE);
40  cols2 = internal::random<int>(2,EIGEN_TEST_MAX_SIZE);
41  }
42  else
43  {
44  cols2 = cols = MatrixType::ColsAtCompileTime;
45  }
46 
47  enum {
48  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
49  ColsAtCompileTime = MatrixType::ColsAtCompileTime
50  };
51  typedef typename internal::kernel_retval_base<FullPivLU<MatrixType> >::ReturnType KernelMatrixType;
52  typedef typename internal::image_retval_base<FullPivLU<MatrixType> >::ReturnType ImageMatrixType;
54  CMatrixType;
56  RMatrixType;
57 
58  Index rank = internal::random<Index>(1, (std::min)(rows, cols)-1);
59 
60  // The image of the zero matrix should consist of a single (zero) column vector
61  VERIFY((MatrixType::Zero(rows,cols).fullPivLu().image(MatrixType::Zero(rows,cols)).cols() == 1));
62 
63  // The kernel of the zero matrix is the entire space, and thus is an invertible matrix of dimensions cols.
64  KernelMatrixType kernel = MatrixType::Zero(rows,cols).fullPivLu().kernel();
65  VERIFY((kernel.fullPivLu().isInvertible()));
66 
67  MatrixType m1(rows, cols), m3(rows, cols2);
68  CMatrixType m2(cols, cols2);
70 
72 
73  // The special value 0.01 below works well in tests. Keep in mind that we're only computing the rank
74  // of singular values are either 0 or 1.
75  // So it's not clear at all that the epsilon should play any role there.
76  lu.setThreshold(RealScalar(0.01));
77  lu.compute(m1);
78 
79  MatrixType u(rows,cols);
80  u = lu.matrixLU().template triangularView<Upper>();
81  RMatrixType l = RMatrixType::Identity(rows,rows);
82  l.block(0,0,rows,(std::min)(rows,cols)).template triangularView<StrictlyLower>()
83  = lu.matrixLU().block(0,0,rows,(std::min)(rows,cols));
84 
85  VERIFY_IS_APPROX(lu.permutationP() * m1 * lu.permutationQ(), l*u);
86 
87  KernelMatrixType m1kernel = lu.kernel();
88  ImageMatrixType m1image = lu.image(m1);
89 
90  VERIFY_IS_APPROX(m1, lu.reconstructedMatrix());
91  VERIFY(rank == lu.rank());
92  VERIFY(cols - lu.rank() == lu.dimensionOfKernel());
93  VERIFY(!lu.isInjective());
94  VERIFY(!lu.isInvertible());
95  VERIFY(!lu.isSurjective());
96  VERIFY_IS_MUCH_SMALLER_THAN((m1 * m1kernel), m1);
97  VERIFY(m1image.fullPivLu().rank() == rank);
98  VERIFY_IS_APPROX(m1 * m1.adjoint() * m1image, m1image);
99 
100  check_solverbase<CMatrixType, MatrixType>(m1, lu, rows, cols, cols2);
101 
102  m2 = CMatrixType::Random(cols,cols2);
103  m3 = m1*m2;
104  m2 = CMatrixType::Random(cols,cols2);
105  // test that the code, which does resize(), may be applied to an xpr
106  m2.block(0,0,m2.rows(),m2.cols()) = lu.solve(m3);
108 }
109 
110 template<typename MatrixType> void lu_invertible()
111 {
112  /* this test covers the following files:
113  FullPivLU.h
114  */
116  Index size = MatrixType::RowsAtCompileTime;
117  if( size==Dynamic)
118  size = internal::random<Index>(1,EIGEN_TEST_MAX_SIZE);
119 
122  lu.setThreshold(RealScalar(0.01));
123  do {
124  m1 = MatrixType::Random(size,size);
125  lu.compute(m1);
126  } while(!lu.isInvertible());
127 
128  VERIFY_IS_APPROX(m1, lu.reconstructedMatrix());
129  VERIFY(0 == lu.dimensionOfKernel());
130  VERIFY(lu.kernel().cols() == 1); // the kernel() should consist of a single (zero) column vector
131  VERIFY(size == lu.rank());
132  VERIFY(lu.isInjective());
133  VERIFY(lu.isSurjective());
134  VERIFY(lu.isInvertible());
135  VERIFY(lu.image(m1).fullPivLu().isInvertible());
136 
137  check_solverbase<MatrixType, MatrixType>(m1, lu, size, size, size);
138 
139  MatrixType m1_inverse = lu.inverse();
140  m3 = MatrixType::Random(size,size);
141  m2 = lu.solve(m3);
142  VERIFY_IS_APPROX(m2, m1_inverse*m3);
143 
144  RealScalar rcond = (RealScalar(1) / matrix_l1_norm(m1)) / matrix_l1_norm(m1_inverse);
145  const RealScalar rcond_est = lu.rcond();
146  // Verify that the estimated condition number is within a factor of 10 of the
147  // truth.
148  VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10);
149 
150  // Regression test for Bug 302
151  MatrixType m4 = MatrixType::Random(size,size);
152  VERIFY_IS_APPROX(lu.solve(m3*m4), lu.solve(m3)*m4);
153 }
154 
155 template<typename MatrixType> void lu_partial_piv(Index size = MatrixType::ColsAtCompileTime)
156 {
157  /* this test covers the following files:
158  PartialPivLU.h
159  */
161 
163  m1.setRandom();
165 
166  STATIC_CHECK(( internal::is_same<typename PartialPivLU<MatrixType>::StorageIndex,int>::value ));
167 
169 
170  check_solverbase<MatrixType, MatrixType>(m1, plu, size, size, size);
171 
172  MatrixType m1_inverse = plu.inverse();
173  m3 = MatrixType::Random(size,size);
174  m2 = plu.solve(m3);
175  VERIFY_IS_APPROX(m2, m1_inverse*m3);
176 
177  RealScalar rcond = (RealScalar(1) / matrix_l1_norm(m1)) / matrix_l1_norm(m1_inverse);
178  const RealScalar rcond_est = plu.rcond();
179  // Verify that the estimate is within a factor of 10 of the truth.
180  VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10);
181 }
182 
183 template<typename MatrixType> void lu_verify_assert()
184 {
185  MatrixType tmp;
186 
188  VERIFY_RAISES_ASSERT(lu.matrixLU())
189  VERIFY_RAISES_ASSERT(lu.permutationP())
190  VERIFY_RAISES_ASSERT(lu.permutationQ())
191  VERIFY_RAISES_ASSERT(lu.kernel())
192  VERIFY_RAISES_ASSERT(lu.image(tmp))
193  VERIFY_RAISES_ASSERT(lu.solve(tmp))
194  VERIFY_RAISES_ASSERT(lu.transpose().solve(tmp))
195  VERIFY_RAISES_ASSERT(lu.adjoint().solve(tmp))
196  VERIFY_RAISES_ASSERT(lu.determinant())
197  VERIFY_RAISES_ASSERT(lu.rank())
198  VERIFY_RAISES_ASSERT(lu.dimensionOfKernel())
199  VERIFY_RAISES_ASSERT(lu.isInjective())
200  VERIFY_RAISES_ASSERT(lu.isSurjective())
201  VERIFY_RAISES_ASSERT(lu.isInvertible())
202  VERIFY_RAISES_ASSERT(lu.inverse())
203 
207  VERIFY_RAISES_ASSERT(plu.solve(tmp))
208  VERIFY_RAISES_ASSERT(plu.transpose().solve(tmp))
209  VERIFY_RAISES_ASSERT(plu.adjoint().solve(tmp))
212 }
213 
215 {
216  for(int i = 0; i < g_repeat; i++) {
217  CALL_SUBTEST_1( lu_non_invertible<Matrix3f>() );
218  CALL_SUBTEST_1( lu_invertible<Matrix3f>() );
219  CALL_SUBTEST_1( lu_verify_assert<Matrix3f>() );
220  CALL_SUBTEST_1( lu_partial_piv<Matrix3f>() );
221 
224  CALL_SUBTEST_2( lu_partial_piv<Matrix2d>() );
225  CALL_SUBTEST_2( lu_partial_piv<Matrix4d>() );
227 
228  CALL_SUBTEST_3( lu_non_invertible<MatrixXf>() );
229  CALL_SUBTEST_3( lu_invertible<MatrixXf>() );
230  CALL_SUBTEST_3( lu_verify_assert<MatrixXf>() );
231 
232  CALL_SUBTEST_4( lu_non_invertible<MatrixXd>() );
233  CALL_SUBTEST_4( lu_invertible<MatrixXd>() );
234  CALL_SUBTEST_4( lu_partial_piv<MatrixXd>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE)) );
235  CALL_SUBTEST_4( lu_verify_assert<MatrixXd>() );
236 
237  CALL_SUBTEST_5( lu_non_invertible<MatrixXcf>() );
238  CALL_SUBTEST_5( lu_invertible<MatrixXcf>() );
239  CALL_SUBTEST_5( lu_verify_assert<MatrixXcf>() );
240 
241  CALL_SUBTEST_6( lu_non_invertible<MatrixXcd>() );
242  CALL_SUBTEST_6( lu_invertible<MatrixXcd>() );
243  CALL_SUBTEST_6( lu_partial_piv<MatrixXcd>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE)) );
244  CALL_SUBTEST_6( lu_verify_assert<MatrixXcd>() );
245 
247 
248  // Test problem size constructors
251  }
252 }
VERIFY_IS_MUCH_SMALLER_THAN
#define VERIFY_IS_MUCH_SMALLER_THAN(a, b)
Definition: main.h:390
Eigen::PartialPivLU
LU decomposition of a matrix with partial pivoting, and related features.
Definition: ForwardDeclarations.h:269
MatrixType
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
m1
Matrix3d m1
Definition: IOFormat.cpp:2
matrix_l1_norm
MatrixType::RealScalar matrix_l1_norm(const MatrixType &m)
Definition: test/lu.cpp:16
CALL_SUBTEST_9
#define CALL_SUBTEST_9(FUNC)
Definition: split_test_helper.h:52
Eigen::FullPivLU
LU decomposition of a matrix with complete pivoting, and related features.
Definition: ForwardDeclarations.h:268
lu_invertible
void lu_invertible()
Definition: test/lu.cpp:110
rows
int rows
Definition: Tutorial_commainit_02.cpp:1
VERIFY_RAISES_ASSERT
#define VERIFY_RAISES_ASSERT(a)
Definition: main.h:340
size
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
lu_verify_assert
void lu_verify_assert()
Definition: test/lu.cpp:183
CALL_SUBTEST_4
#define CALL_SUBTEST_4(FUNC)
Definition: split_test_helper.h:22
Eigen::PartialPivLU::determinant
Scalar determinant() const
Definition: PartialPivLU.h:552
m2
MatrixType m2(n_dims)
CALL_SUBTEST_3
#define CALL_SUBTEST_3(FUNC)
Definition: split_test_helper.h:16
CALL_SUBTEST_1
#define CALL_SUBTEST_1(FUNC)
Definition: split_test_helper.h:4
Eigen::createRandomPIMatrixOfRank
void createRandomPIMatrixOfRank(Index desired_rank, Index rows, Index cols, MatrixType &m)
Definition: main.h:653
Eigen::Dynamic
const int Dynamic
Definition: Constants.h:22
l
static const Line3 l(Rot3(), 1, 1)
Eigen::PartialPivLU::matrixLU
const MatrixType & matrixLU() const
Definition: PartialPivLU.h:143
CALL_SUBTEST_5
#define CALL_SUBTEST_5(FUNC)
Definition: split_test_helper.h:28
Eigen::g_repeat
static int g_repeat
Definition: main.h:169
m
Matrix3f m
Definition: AngleAxis_mimic_euler.cpp:1
lu_non_invertible
void lu_non_invertible()
Definition: test/lu.cpp:20
CALL_SUBTEST_6
#define CALL_SUBTEST_6(FUNC)
Definition: split_test_helper.h:34
CALL_SUBTEST_2
#define CALL_SUBTEST_2(FUNC)
Definition: split_test_helper.h:10
lu_partial_piv
void lu_partial_piv(Index size=MatrixType::ColsAtCompileTime)
Definition: test/lu.cpp:155
VERIFY_IS_APPROX
#define VERIFY_IS_APPROX(a, b)
Definition: integer_types.cpp:15
RealScalar
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:47
lu
cout<< "Here is the matrix m:"<< endl<< m<< endl;Eigen::FullPivLU< Matrix5x3 > lu(m)
m3
static const DiscreteKey m3(M(3), 2)
solverbase.h
Eigen::PartialPivLU::rcond
RealScalar rcond() const
Definition: PartialPivLU.h:183
Eigen::PartialPivLU::inverse
const Inverse< PartialPivLU > inverse() const
Definition: PartialPivLU.h:196
main.h
std
Definition: BFloat16.h:88
EIGEN_TEST_MAX_SIZE
#define EIGEN_TEST_MAX_SIZE
Definition: boostmultiprec.cpp:16
Eigen::PartialPivLU::reconstructedMatrix
MatrixType reconstructedMatrix() const
Definition: PartialPivLU.h:562
min
#define min(a, b)
Definition: datatypes.h:19
Eigen::Matrix
The matrix class, also used for vectors and row-vectors.
Definition: 3rdparty/Eigen/Eigen/src/Core/Matrix.h:178
EIGEN_DECLARE_TEST
EIGEN_DECLARE_TEST(lu)
Definition: test/lu.cpp:214
cols
int cols
Definition: Tutorial_commainit_02.cpp:1
CALL_SUBTEST_7
#define CALL_SUBTEST_7(FUNC)
Definition: split_test_helper.h:40
Eigen::NumTraits
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:232
test_callbacks.value
value
Definition: test_callbacks.py:160
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9
STATIC_CHECK
#define STATIC_CHECK(COND)
Definition: main.h:397
Eigen::PartialPivLU::permutationP
const PermutationType & permutationP() const
Definition: PartialPivLU.h:151
VERIFY
#define VERIFY(a)
Definition: main.h:380
Eigen::Index
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74


gtsam
Author(s):
autogenerated on Fri Nov 1 2024 03:33:10