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
lu_non_invertible()14 template<typename MatrixType> void lu_non_invertible()
15 {
16 typedef typename MatrixType::Index Index;
17 typedef typename MatrixType::RealScalar RealScalar;
18 /* this test covers the following files:
19 LU.h
20 */
21 Index rows, cols, cols2;
22 if(MatrixType::RowsAtCompileTime==Dynamic)
23 {
24 rows = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE);
25 }
26 else
27 {
28 rows = MatrixType::RowsAtCompileTime;
29 }
30 if(MatrixType::ColsAtCompileTime==Dynamic)
31 {
32 cols = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE);
33 cols2 = internal::random<int>(2,EIGEN_TEST_MAX_SIZE);
34 }
35 else
36 {
37 cols2 = cols = MatrixType::ColsAtCompileTime;
38 }
39
40 enum {
41 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
42 ColsAtCompileTime = MatrixType::ColsAtCompileTime
43 };
44 typedef typename internal::kernel_retval_base<FullPivLU<MatrixType> >::ReturnType KernelMatrixType;
45 typedef typename internal::image_retval_base<FullPivLU<MatrixType> >::ReturnType ImageMatrixType;
46 typedef Matrix<typename MatrixType::Scalar, ColsAtCompileTime, ColsAtCompileTime>
47 CMatrixType;
48 typedef Matrix<typename MatrixType::Scalar, RowsAtCompileTime, RowsAtCompileTime>
49 RMatrixType;
50
51 Index rank = internal::random<Index>(1, (std::min)(rows, cols)-1);
52
53 // The image of the zero matrix should consist of a single (zero) column vector
54 VERIFY((MatrixType::Zero(rows,cols).fullPivLu().image(MatrixType::Zero(rows,cols)).cols() == 1));
55
56 MatrixType m1(rows, cols), m3(rows, cols2);
57 CMatrixType m2(cols, cols2);
58 createRandomPIMatrixOfRank(rank, rows, cols, m1);
59
60 FullPivLU<MatrixType> lu;
61
62 // The special value 0.01 below works well in tests. Keep in mind that we're only computing the rank
63 // of singular values are either 0 or 1.
64 // So it's not clear at all that the epsilon should play any role there.
65 lu.setThreshold(RealScalar(0.01));
66 lu.compute(m1);
67
68 MatrixType u(rows,cols);
69 u = lu.matrixLU().template triangularView<Upper>();
70 RMatrixType l = RMatrixType::Identity(rows,rows);
71 l.block(0,0,rows,(std::min)(rows,cols)).template triangularView<StrictlyLower>()
72 = lu.matrixLU().block(0,0,rows,(std::min)(rows,cols));
73
74 VERIFY_IS_APPROX(lu.permutationP() * m1 * lu.permutationQ(), l*u);
75
76 KernelMatrixType m1kernel = lu.kernel();
77 ImageMatrixType m1image = lu.image(m1);
78
79 VERIFY_IS_APPROX(m1, lu.reconstructedMatrix());
80 VERIFY(rank == lu.rank());
81 VERIFY(cols - lu.rank() == lu.dimensionOfKernel());
82 VERIFY(!lu.isInjective());
83 VERIFY(!lu.isInvertible());
84 VERIFY(!lu.isSurjective());
85 VERIFY((m1 * m1kernel).isMuchSmallerThan(m1));
86 VERIFY(m1image.fullPivLu().rank() == rank);
87 VERIFY_IS_APPROX(m1 * m1.adjoint() * m1image, m1image);
88
89 m2 = CMatrixType::Random(cols,cols2);
90 m3 = m1*m2;
91 m2 = CMatrixType::Random(cols,cols2);
92 // test that the code, which does resize(), may be applied to an xpr
93 m2.block(0,0,m2.rows(),m2.cols()) = lu.solve(m3);
94 VERIFY_IS_APPROX(m3, m1*m2);
95 }
96
lu_invertible()97 template<typename MatrixType> void lu_invertible()
98 {
99 /* this test covers the following files:
100 LU.h
101 */
102 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
103 DenseIndex size = MatrixType::RowsAtCompileTime;
104 if( size==Dynamic)
105 size = internal::random<DenseIndex>(1,EIGEN_TEST_MAX_SIZE);
106
107 MatrixType m1(size, size), m2(size, size), m3(size, size);
108 FullPivLU<MatrixType> lu;
109 lu.setThreshold(RealScalar(0.01));
110 do {
111 m1 = MatrixType::Random(size,size);
112 lu.compute(m1);
113 } while(!lu.isInvertible());
114
115 VERIFY_IS_APPROX(m1, lu.reconstructedMatrix());
116 VERIFY(0 == lu.dimensionOfKernel());
117 VERIFY(lu.kernel().cols() == 1); // the kernel() should consist of a single (zero) column vector
118 VERIFY(size == lu.rank());
119 VERIFY(lu.isInjective());
120 VERIFY(lu.isSurjective());
121 VERIFY(lu.isInvertible());
122 VERIFY(lu.image(m1).fullPivLu().isInvertible());
123 m3 = MatrixType::Random(size,size);
124 m2 = lu.solve(m3);
125 VERIFY_IS_APPROX(m3, m1*m2);
126 VERIFY_IS_APPROX(m2, lu.inverse()*m3);
127
128 // Regression test for Bug 302
129 MatrixType m4 = MatrixType::Random(size,size);
130 VERIFY_IS_APPROX(lu.solve(m3*m4), lu.solve(m3)*m4);
131 }
132
lu_partial_piv()133 template<typename MatrixType> void lu_partial_piv()
134 {
135 /* this test covers the following files:
136 PartialPivLU.h
137 */
138 typedef typename MatrixType::Index Index;
139 Index rows = internal::random<Index>(1,4);
140 Index cols = rows;
141
142 MatrixType m1(cols, rows);
143 m1.setRandom();
144 PartialPivLU<MatrixType> plu(m1);
145
146 VERIFY_IS_APPROX(m1, plu.reconstructedMatrix());
147 }
148
lu_verify_assert()149 template<typename MatrixType> void lu_verify_assert()
150 {
151 MatrixType tmp;
152
153 FullPivLU<MatrixType> lu;
154 VERIFY_RAISES_ASSERT(lu.matrixLU())
155 VERIFY_RAISES_ASSERT(lu.permutationP())
156 VERIFY_RAISES_ASSERT(lu.permutationQ())
157 VERIFY_RAISES_ASSERT(lu.kernel())
158 VERIFY_RAISES_ASSERT(lu.image(tmp))
159 VERIFY_RAISES_ASSERT(lu.solve(tmp))
160 VERIFY_RAISES_ASSERT(lu.determinant())
161 VERIFY_RAISES_ASSERT(lu.rank())
162 VERIFY_RAISES_ASSERT(lu.dimensionOfKernel())
163 VERIFY_RAISES_ASSERT(lu.isInjective())
164 VERIFY_RAISES_ASSERT(lu.isSurjective())
165 VERIFY_RAISES_ASSERT(lu.isInvertible())
166 VERIFY_RAISES_ASSERT(lu.inverse())
167
168 PartialPivLU<MatrixType> plu;
169 VERIFY_RAISES_ASSERT(plu.matrixLU())
170 VERIFY_RAISES_ASSERT(plu.permutationP())
171 VERIFY_RAISES_ASSERT(plu.solve(tmp))
172 VERIFY_RAISES_ASSERT(plu.determinant())
173 VERIFY_RAISES_ASSERT(plu.inverse())
174 }
175
test_lu()176 void test_lu()
177 {
178 for(int i = 0; i < g_repeat; i++) {
179 CALL_SUBTEST_1( lu_non_invertible<Matrix3f>() );
180 CALL_SUBTEST_1( lu_invertible<Matrix3f>() );
181 CALL_SUBTEST_1( lu_verify_assert<Matrix3f>() );
182
183 CALL_SUBTEST_2( (lu_non_invertible<Matrix<double, 4, 6> >()) );
184 CALL_SUBTEST_2( (lu_verify_assert<Matrix<double, 4, 6> >()) );
185
186 CALL_SUBTEST_3( lu_non_invertible<MatrixXf>() );
187 CALL_SUBTEST_3( lu_invertible<MatrixXf>() );
188 CALL_SUBTEST_3( lu_verify_assert<MatrixXf>() );
189
190 CALL_SUBTEST_4( lu_non_invertible<MatrixXd>() );
191 CALL_SUBTEST_4( lu_invertible<MatrixXd>() );
192 CALL_SUBTEST_4( lu_partial_piv<MatrixXd>() );
193 CALL_SUBTEST_4( lu_verify_assert<MatrixXd>() );
194
195 CALL_SUBTEST_5( lu_non_invertible<MatrixXcf>() );
196 CALL_SUBTEST_5( lu_invertible<MatrixXcf>() );
197 CALL_SUBTEST_5( lu_verify_assert<MatrixXcf>() );
198
199 CALL_SUBTEST_6( lu_non_invertible<MatrixXcd>() );
200 CALL_SUBTEST_6( lu_invertible<MatrixXcd>() );
201 CALL_SUBTEST_6( lu_partial_piv<MatrixXcd>() );
202 CALL_SUBTEST_6( lu_verify_assert<MatrixXcd>() );
203
204 CALL_SUBTEST_7(( lu_non_invertible<Matrix<float,Dynamic,16> >() ));
205
206 // Test problem size constructors
207 CALL_SUBTEST_9( PartialPivLU<MatrixXf>(10) );
208 CALL_SUBTEST_9( FullPivLU<MatrixXf>(10, 20); );
209 }
210 }
211