Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 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 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 14 | template<typename MatrixType> |
| 15 | typename MatrixType::RealScalar matrix_l1_norm(const MatrixType& m) { |
| 16 | return m.cwiseAbs().colwise().sum().maxCoeff(); |
| 17 | } |
| 18 | |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 19 | template<typename MatrixType> void lu_non_invertible() |
| 20 | { |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 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; |
| 50 | typedef Matrix<typename MatrixType::Scalar, ColsAtCompileTime, ColsAtCompileTime> |
| 51 | CMatrixType; |
| 52 | typedef Matrix<typename MatrixType::Scalar, RowsAtCompileTime, RowsAtCompileTime> |
| 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 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 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 | |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 64 | MatrixType m1(rows, cols), m3(rows, cols2); |
| 65 | CMatrixType m2(cols, cols2); |
| 66 | createRandomPIMatrixOfRank(rank, rows, cols, m1); |
| 67 | |
| 68 | FullPivLU<MatrixType> lu; |
| 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 | |
| 87 | VERIFY_IS_APPROX(m1, lu.reconstructedMatrix()); |
| 88 | VERIFY(rank == lu.rank()); |
| 89 | VERIFY(cols - lu.rank() == lu.dimensionOfKernel()); |
| 90 | VERIFY(!lu.isInjective()); |
| 91 | VERIFY(!lu.isInvertible()); |
| 92 | VERIFY(!lu.isSurjective()); |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 93 | VERIFY_IS_MUCH_SMALLER_THAN((m1 * m1kernel), m1); |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 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); |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 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); |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 123 | } |
| 124 | |
| 125 | template<typename MatrixType> void lu_invertible() |
| 126 | { |
| 127 | /* this test covers the following files: |
| 128 | LU.h |
| 129 | */ |
| 130 | typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar; |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 131 | Index size = MatrixType::RowsAtCompileTime; |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 132 | if( size==Dynamic) |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 133 | size = internal::random<Index>(1,EIGEN_TEST_MAX_SIZE); |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 134 | |
| 135 | MatrixType m1(size, size), m2(size, size), m3(size, size); |
| 136 | FullPivLU<MatrixType> lu; |
| 137 | lu.setThreshold(RealScalar(0.01)); |
| 138 | do { |
| 139 | m1 = MatrixType::Random(size,size); |
| 140 | lu.compute(m1); |
| 141 | } while(!lu.isInvertible()); |
| 142 | |
| 143 | VERIFY_IS_APPROX(m1, lu.reconstructedMatrix()); |
| 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); |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 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); |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 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 | */ |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 187 | typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar; |
| 188 | Index size = internal::random<Index>(1,4); |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 189 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 190 | MatrixType m1(size, size), m2(size, size), m3(size, size); |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 191 | m1.setRandom(); |
| 192 | PartialPivLU<MatrixType> plu(m1); |
| 193 | |
| 194 | VERIFY_IS_APPROX(m1, plu.reconstructedMatrix()); |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 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); |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 220 | } |
| 221 | |
| 222 | template<typename MatrixType> void lu_verify_assert() |
| 223 | { |
| 224 | MatrixType tmp; |
| 225 | |
| 226 | FullPivLU<MatrixType> lu; |
| 227 | VERIFY_RAISES_ASSERT(lu.matrixLU()) |
| 228 | VERIFY_RAISES_ASSERT(lu.permutationP()) |
| 229 | VERIFY_RAISES_ASSERT(lu.permutationQ()) |
| 230 | VERIFY_RAISES_ASSERT(lu.kernel()) |
| 231 | VERIFY_RAISES_ASSERT(lu.image(tmp)) |
| 232 | VERIFY_RAISES_ASSERT(lu.solve(tmp)) |
| 233 | VERIFY_RAISES_ASSERT(lu.determinant()) |
| 234 | VERIFY_RAISES_ASSERT(lu.rank()) |
| 235 | VERIFY_RAISES_ASSERT(lu.dimensionOfKernel()) |
| 236 | VERIFY_RAISES_ASSERT(lu.isInjective()) |
| 237 | VERIFY_RAISES_ASSERT(lu.isSurjective()) |
| 238 | VERIFY_RAISES_ASSERT(lu.isInvertible()) |
| 239 | VERIFY_RAISES_ASSERT(lu.inverse()) |
| 240 | |
| 241 | PartialPivLU<MatrixType> plu; |
| 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 | } |