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