permutationmatrices.cpp 6.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 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. #define TEST_ENABLE_TEMPORARY_TRACKING
  10. #include "main.h"
  11. using namespace std;
  12. template<typename MatrixType> void permutationmatrices(const MatrixType& m)
  13. {
  14. typedef typename MatrixType::Scalar Scalar;
  15. enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime,
  16. Options = MatrixType::Options };
  17. typedef PermutationMatrix<Rows> LeftPermutationType;
  18. typedef Transpositions<Rows> LeftTranspositionsType;
  19. typedef Matrix<int, Rows, 1> LeftPermutationVectorType;
  20. typedef Map<LeftPermutationType> MapLeftPerm;
  21. typedef PermutationMatrix<Cols> RightPermutationType;
  22. typedef Transpositions<Cols> RightTranspositionsType;
  23. typedef Matrix<int, Cols, 1> RightPermutationVectorType;
  24. typedef Map<RightPermutationType> MapRightPerm;
  25. Index rows = m.rows();
  26. Index cols = m.cols();
  27. MatrixType m_original = MatrixType::Random(rows,cols);
  28. LeftPermutationVectorType lv;
  29. randomPermutationVector(lv, rows);
  30. LeftPermutationType lp(lv);
  31. RightPermutationVectorType rv;
  32. randomPermutationVector(rv, cols);
  33. RightPermutationType rp(rv);
  34. LeftTranspositionsType lt(lv);
  35. RightTranspositionsType rt(rv);
  36. MatrixType m_permuted = MatrixType::Random(rows,cols);
  37. VERIFY_EVALUATION_COUNT(m_permuted = lp * m_original * rp, 1); // 1 temp for sub expression "lp * m_original"
  38. for (int i=0; i<rows; i++)
  39. for (int j=0; j<cols; j++)
  40. VERIFY_IS_APPROX(m_permuted(lv(i),j), m_original(i,rv(j)));
  41. Matrix<Scalar,Rows,Rows> lm(lp);
  42. Matrix<Scalar,Cols,Cols> rm(rp);
  43. VERIFY_IS_APPROX(m_permuted, lm*m_original*rm);
  44. m_permuted = m_original;
  45. VERIFY_EVALUATION_COUNT(m_permuted = lp * m_permuted * rp, 1);
  46. VERIFY_IS_APPROX(m_permuted, lm*m_original*rm);
  47. LeftPermutationType lpi;
  48. lpi = lp.inverse();
  49. VERIFY_IS_APPROX(lpi*m_permuted,lp.inverse()*m_permuted);
  50. VERIFY_IS_APPROX(lp.inverse()*m_permuted*rp.inverse(), m_original);
  51. VERIFY_IS_APPROX(lv.asPermutation().inverse()*m_permuted*rv.asPermutation().inverse(), m_original);
  52. VERIFY_IS_APPROX(MapLeftPerm(lv.data(),lv.size()).inverse()*m_permuted*MapRightPerm(rv.data(),rv.size()).inverse(), m_original);
  53. VERIFY((lp*lp.inverse()).toDenseMatrix().isIdentity());
  54. VERIFY((lv.asPermutation()*lv.asPermutation().inverse()).toDenseMatrix().isIdentity());
  55. VERIFY((MapLeftPerm(lv.data(),lv.size())*MapLeftPerm(lv.data(),lv.size()).inverse()).toDenseMatrix().isIdentity());
  56. LeftPermutationVectorType lv2;
  57. randomPermutationVector(lv2, rows);
  58. LeftPermutationType lp2(lv2);
  59. Matrix<Scalar,Rows,Rows> lm2(lp2);
  60. VERIFY_IS_APPROX((lp*lp2).toDenseMatrix().template cast<Scalar>(), lm*lm2);
  61. VERIFY_IS_APPROX((lv.asPermutation()*lv2.asPermutation()).toDenseMatrix().template cast<Scalar>(), lm*lm2);
  62. VERIFY_IS_APPROX((MapLeftPerm(lv.data(),lv.size())*MapLeftPerm(lv2.data(),lv2.size())).toDenseMatrix().template cast<Scalar>(), lm*lm2);
  63. LeftPermutationType identityp;
  64. identityp.setIdentity(rows);
  65. VERIFY_IS_APPROX(m_original, identityp*m_original);
  66. // check inplace permutations
  67. m_permuted = m_original;
  68. VERIFY_EVALUATION_COUNT(m_permuted.noalias()= lp.inverse() * m_permuted, 1); // 1 temp to allocate the mask
  69. VERIFY_IS_APPROX(m_permuted, lp.inverse()*m_original);
  70. m_permuted = m_original;
  71. VERIFY_EVALUATION_COUNT(m_permuted.noalias() = m_permuted * rp.inverse(), 1); // 1 temp to allocate the mask
  72. VERIFY_IS_APPROX(m_permuted, m_original*rp.inverse());
  73. m_permuted = m_original;
  74. VERIFY_EVALUATION_COUNT(m_permuted.noalias() = lp * m_permuted, 1); // 1 temp to allocate the mask
  75. VERIFY_IS_APPROX(m_permuted, lp*m_original);
  76. m_permuted = m_original;
  77. VERIFY_EVALUATION_COUNT(m_permuted.noalias() = m_permuted * rp, 1); // 1 temp to allocate the mask
  78. VERIFY_IS_APPROX(m_permuted, m_original*rp);
  79. if(rows>1 && cols>1)
  80. {
  81. lp2 = lp;
  82. Index i = internal::random<Index>(0, rows-1);
  83. Index j;
  84. do j = internal::random<Index>(0, rows-1); while(j==i);
  85. lp2.applyTranspositionOnTheLeft(i, j);
  86. lm = lp;
  87. lm.row(i).swap(lm.row(j));
  88. VERIFY_IS_APPROX(lm, lp2.toDenseMatrix().template cast<Scalar>());
  89. RightPermutationType rp2 = rp;
  90. i = internal::random<Index>(0, cols-1);
  91. do j = internal::random<Index>(0, cols-1); while(j==i);
  92. rp2.applyTranspositionOnTheRight(i, j);
  93. rm = rp;
  94. rm.col(i).swap(rm.col(j));
  95. VERIFY_IS_APPROX(rm, rp2.toDenseMatrix().template cast<Scalar>());
  96. }
  97. {
  98. // simple compilation check
  99. Matrix<Scalar, Cols, Cols> A = rp;
  100. Matrix<Scalar, Cols, Cols> B = rp.transpose();
  101. VERIFY_IS_APPROX(A, B.transpose());
  102. }
  103. m_permuted = m_original;
  104. lp = lt;
  105. rp = rt;
  106. VERIFY_EVALUATION_COUNT(m_permuted = lt * m_permuted * rt, 1);
  107. VERIFY_IS_APPROX(m_permuted, lp*m_original*rp.transpose());
  108. VERIFY_IS_APPROX(lt.inverse()*m_permuted*rt.inverse(), m_original);
  109. // Check inplace transpositions
  110. m_permuted = m_original;
  111. VERIFY_IS_APPROX(m_permuted = lt * m_permuted, lp * m_original);
  112. m_permuted = m_original;
  113. VERIFY_IS_APPROX(m_permuted = lt.inverse() * m_permuted, lp.inverse() * m_original);
  114. m_permuted = m_original;
  115. VERIFY_IS_APPROX(m_permuted = m_permuted * rt, m_original * rt);
  116. m_permuted = m_original;
  117. VERIFY_IS_APPROX(m_permuted = m_permuted * rt.inverse(), m_original * rt.inverse());
  118. }
  119. template<typename T>
  120. void bug890()
  121. {
  122. typedef Matrix<T, Dynamic, Dynamic> MatrixType;
  123. typedef Matrix<T, Dynamic, 1> VectorType;
  124. typedef Stride<Dynamic,Dynamic> S;
  125. typedef Map<MatrixType, Aligned, S> MapType;
  126. typedef PermutationMatrix<Dynamic> Perm;
  127. VectorType v1(2), v2(2), op(4), rhs(2);
  128. v1 << 666,667;
  129. op << 1,0,0,1;
  130. rhs << 42,42;
  131. Perm P(2);
  132. P.indices() << 1, 0;
  133. MapType(v1.data(),2,1,S(1,1)) = P * MapType(rhs.data(),2,1,S(1,1));
  134. VERIFY_IS_APPROX(v1, (P * rhs).eval());
  135. MapType(v1.data(),2,1,S(1,1)) = P.inverse() * MapType(rhs.data(),2,1,S(1,1));
  136. VERIFY_IS_APPROX(v1, (P.inverse() * rhs).eval());
  137. }
  138. EIGEN_DECLARE_TEST(permutationmatrices)
  139. {
  140. for(int i = 0; i < g_repeat; i++) {
  141. CALL_SUBTEST_1( permutationmatrices(Matrix<float, 1, 1>()) );
  142. CALL_SUBTEST_2( permutationmatrices(Matrix3f()) );
  143. CALL_SUBTEST_3( permutationmatrices(Matrix<double,3,3,RowMajor>()) );
  144. CALL_SUBTEST_4( permutationmatrices(Matrix4d()) );
  145. CALL_SUBTEST_5( permutationmatrices(Matrix<double,40,60>()) );
  146. CALL_SUBTEST_6( permutationmatrices(Matrix<double,Dynamic,Dynamic,RowMajor>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  147. CALL_SUBTEST_7( permutationmatrices(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  148. }
  149. CALL_SUBTEST_5( bug890<double>() );
  150. }