inverse.cpp 4.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. // Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
  6. //
  7. // This Source Code Form is subject to the terms of the Mozilla
  8. // Public License v. 2.0. If a copy of the MPL was not distributed
  9. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  10. #include "main.h"
  11. #include <Eigen/LU>
  12. template<typename MatrixType>
  13. void inverse_for_fixed_size(const MatrixType&, typename internal::enable_if<MatrixType::SizeAtCompileTime==Dynamic>::type* = 0)
  14. {
  15. }
  16. template<typename MatrixType>
  17. void inverse_for_fixed_size(const MatrixType& m1, typename internal::enable_if<MatrixType::SizeAtCompileTime!=Dynamic>::type* = 0)
  18. {
  19. using std::abs;
  20. MatrixType m2, identity = MatrixType::Identity();
  21. typedef typename MatrixType::Scalar Scalar;
  22. typedef typename NumTraits<Scalar>::Real RealScalar;
  23. typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
  24. //computeInverseAndDetWithCheck tests
  25. //First: an invertible matrix
  26. bool invertible;
  27. Scalar det;
  28. m2.setZero();
  29. m1.computeInverseAndDetWithCheck(m2, det, invertible);
  30. VERIFY(invertible);
  31. VERIFY_IS_APPROX(identity, m1*m2);
  32. VERIFY_IS_APPROX(det, m1.determinant());
  33. m2.setZero();
  34. m1.computeInverseWithCheck(m2, invertible);
  35. VERIFY(invertible);
  36. VERIFY_IS_APPROX(identity, m1*m2);
  37. //Second: a rank one matrix (not invertible, except for 1x1 matrices)
  38. VectorType v3 = VectorType::Random();
  39. MatrixType m3 = v3*v3.transpose(), m4;
  40. m3.computeInverseAndDetWithCheck(m4, det, invertible);
  41. VERIFY( m1.rows()==1 ? invertible : !invertible );
  42. VERIFY_IS_MUCH_SMALLER_THAN(abs(det-m3.determinant()), RealScalar(1));
  43. m3.computeInverseWithCheck(m4, invertible);
  44. VERIFY( m1.rows()==1 ? invertible : !invertible );
  45. // check with submatrices
  46. {
  47. Matrix<Scalar, MatrixType::RowsAtCompileTime+1, MatrixType::RowsAtCompileTime+1, MatrixType::Options> m5;
  48. m5.setRandom();
  49. m5.topLeftCorner(m1.rows(),m1.rows()) = m1;
  50. m2 = m5.template topLeftCorner<MatrixType::RowsAtCompileTime,MatrixType::ColsAtCompileTime>().inverse();
  51. VERIFY_IS_APPROX( (m5.template topLeftCorner<MatrixType::RowsAtCompileTime,MatrixType::ColsAtCompileTime>()), m2.inverse() );
  52. }
  53. }
  54. template<typename MatrixType> void inverse(const MatrixType& m)
  55. {
  56. /* this test covers the following files:
  57. Inverse.h
  58. */
  59. Index rows = m.rows();
  60. Index cols = m.cols();
  61. typedef typename MatrixType::Scalar Scalar;
  62. MatrixType m1(rows, cols),
  63. m2(rows, cols),
  64. identity = MatrixType::Identity(rows, rows);
  65. createRandomPIMatrixOfRank(rows,rows,rows,m1);
  66. m2 = m1.inverse();
  67. VERIFY_IS_APPROX(m1, m2.inverse() );
  68. VERIFY_IS_APPROX((Scalar(2)*m2).inverse(), m2.inverse()*Scalar(0.5));
  69. VERIFY_IS_APPROX(identity, m1.inverse() * m1 );
  70. VERIFY_IS_APPROX(identity, m1 * m1.inverse() );
  71. VERIFY_IS_APPROX(m1, m1.inverse().inverse() );
  72. // since for the general case we implement separately row-major and col-major, test that
  73. VERIFY_IS_APPROX(MatrixType(m1.transpose().inverse()), MatrixType(m1.inverse().transpose()));
  74. inverse_for_fixed_size(m1);
  75. // check in-place inversion
  76. if(MatrixType::RowsAtCompileTime>=2 && MatrixType::RowsAtCompileTime<=4)
  77. {
  78. // in-place is forbidden
  79. VERIFY_RAISES_ASSERT(m1 = m1.inverse());
  80. }
  81. else
  82. {
  83. m2 = m1.inverse();
  84. m1 = m1.inverse();
  85. VERIFY_IS_APPROX(m1,m2);
  86. }
  87. }
  88. template<typename Scalar>
  89. void inverse_zerosized()
  90. {
  91. Matrix<Scalar,Dynamic,Dynamic> A(0,0);
  92. {
  93. Matrix<Scalar,0,1> b, x;
  94. x = A.inverse() * b;
  95. }
  96. {
  97. Matrix<Scalar,Dynamic,Dynamic> b(0,1), x;
  98. x = A.inverse() * b;
  99. VERIFY_IS_EQUAL(x.rows(), 0);
  100. VERIFY_IS_EQUAL(x.cols(), 1);
  101. }
  102. }
  103. EIGEN_DECLARE_TEST(inverse)
  104. {
  105. int s = 0;
  106. for(int i = 0; i < g_repeat; i++) {
  107. CALL_SUBTEST_1( inverse(Matrix<double,1,1>()) );
  108. CALL_SUBTEST_2( inverse(Matrix2d()) );
  109. CALL_SUBTEST_3( inverse(Matrix3f()) );
  110. CALL_SUBTEST_4( inverse(Matrix4f()) );
  111. CALL_SUBTEST_4( inverse(Matrix<float,4,4,DontAlign>()) );
  112. s = internal::random<int>(50,320);
  113. CALL_SUBTEST_5( inverse(MatrixXf(s,s)) );
  114. TEST_SET_BUT_UNUSED_VARIABLE(s)
  115. CALL_SUBTEST_5( inverse_zerosized<float>() );
  116. CALL_SUBTEST_5( inverse(MatrixXf(0, 0)) );
  117. CALL_SUBTEST_5( inverse(MatrixXf(1, 1)) );
  118. s = internal::random<int>(25,100);
  119. CALL_SUBTEST_6( inverse(MatrixXcd(s,s)) );
  120. TEST_SET_BUT_UNUSED_VARIABLE(s)
  121. CALL_SUBTEST_7( inverse(Matrix4d()) );
  122. CALL_SUBTEST_7( inverse(Matrix<double,4,4,DontAlign>()) );
  123. CALL_SUBTEST_8( inverse(Matrix4cd()) );
  124. }
  125. }