• 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 // Copyright (C) 2009-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
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 
11 #ifndef EIGEN_PERMUTATIONMATRIX_H
12 #define EIGEN_PERMUTATIONMATRIX_H
13 
14 namespace Eigen {
15 
16 template<int RowCol,typename IndicesType,typename MatrixType, typename StorageKind> class PermutedImpl;
17 
18 /** \class PermutationBase
19   * \ingroup Core_Module
20   *
21   * \brief Base class for permutations
22   *
23   * \param Derived the derived class
24   *
25   * This class is the base class for all expressions representing a permutation matrix,
26   * internally stored as a vector of integers.
27   * The convention followed here is that if \f$ \sigma \f$ is a permutation, the corresponding permutation matrix
28   * \f$ P_\sigma \f$ is such that if \f$ (e_1,\ldots,e_p) \f$ is the canonical basis, we have:
29   *  \f[ P_\sigma(e_i) = e_{\sigma(i)}. \f]
30   * This convention ensures that for any two permutations \f$ \sigma, \tau \f$, we have:
31   *  \f[ P_{\sigma\circ\tau} = P_\sigma P_\tau. \f]
32   *
33   * Permutation matrices are square and invertible.
34   *
35   * Notice that in addition to the member functions and operators listed here, there also are non-member
36   * operator* to multiply any kind of permutation object with any kind of matrix expression (MatrixBase)
37   * on either side.
38   *
39   * \sa class PermutationMatrix, class PermutationWrapper
40   */
41 
42 namespace internal {
43 
44 template<typename PermutationType, typename MatrixType, int Side, bool Transposed=false>
45 struct permut_matrix_product_retval;
46 template<typename PermutationType, typename MatrixType, int Side, bool Transposed=false>
47 struct permut_sparsematrix_product_retval;
48 enum PermPermProduct_t {PermPermProduct};
49 
50 } // end namespace internal
51 
52 template<typename Derived>
53 class PermutationBase : public EigenBase<Derived>
54 {
55     typedef internal::traits<Derived> Traits;
56     typedef EigenBase<Derived> Base;
57   public:
58 
59     #ifndef EIGEN_PARSED_BY_DOXYGEN
60     typedef typename Traits::IndicesType IndicesType;
61     enum {
62       Flags = Traits::Flags,
63       CoeffReadCost = Traits::CoeffReadCost,
64       RowsAtCompileTime = Traits::RowsAtCompileTime,
65       ColsAtCompileTime = Traits::ColsAtCompileTime,
66       MaxRowsAtCompileTime = Traits::MaxRowsAtCompileTime,
67       MaxColsAtCompileTime = Traits::MaxColsAtCompileTime
68     };
69     typedef typename Traits::Scalar Scalar;
70     typedef typename Traits::Index Index;
71     typedef Matrix<Scalar,RowsAtCompileTime,ColsAtCompileTime,0,MaxRowsAtCompileTime,MaxColsAtCompileTime>
72             DenseMatrixType;
73     typedef PermutationMatrix<IndicesType::SizeAtCompileTime,IndicesType::MaxSizeAtCompileTime,Index>
74             PlainPermutationType;
75     using Base::derived;
76     #endif
77 
78     /** Copies the other permutation into *this */
79     template<typename OtherDerived>
80     Derived& operator=(const PermutationBase<OtherDerived>& other)
81     {
82       indices() = other.indices();
83       return derived();
84     }
85 
86     /** Assignment from the Transpositions \a tr */
87     template<typename OtherDerived>
88     Derived& operator=(const TranspositionsBase<OtherDerived>& tr)
89     {
90       setIdentity(tr.size());
91       for(Index k=size()-1; k>=0; --k)
92         applyTranspositionOnTheRight(k,tr.coeff(k));
93       return derived();
94     }
95 
96     #ifndef EIGEN_PARSED_BY_DOXYGEN
97     /** This is a special case of the templated operator=. Its purpose is to
98       * prevent a default operator= from hiding the templated operator=.
99       */
100     Derived& operator=(const PermutationBase& other)
101     {
102       indices() = other.indices();
103       return derived();
104     }
105     #endif
106 
107     /** \returns the number of rows */
rows()108     inline Index rows() const { return Index(indices().size()); }
109 
110     /** \returns the number of columns */
cols()111     inline Index cols() const { return Index(indices().size()); }
112 
113     /** \returns the size of a side of the respective square matrix, i.e., the number of indices */
size()114     inline Index size() const { return Index(indices().size()); }
115 
116     #ifndef EIGEN_PARSED_BY_DOXYGEN
117     template<typename DenseDerived>
evalTo(MatrixBase<DenseDerived> & other)118     void evalTo(MatrixBase<DenseDerived>& other) const
119     {
120       other.setZero();
121       for (int i=0; i<rows();++i)
122         other.coeffRef(indices().coeff(i),i) = typename DenseDerived::Scalar(1);
123     }
124     #endif
125 
126     /** \returns a Matrix object initialized from this permutation matrix. Notice that it
127       * is inefficient to return this Matrix object by value. For efficiency, favor using
128       * the Matrix constructor taking EigenBase objects.
129       */
toDenseMatrix()130     DenseMatrixType toDenseMatrix() const
131     {
132       return derived();
133     }
134 
135     /** const version of indices(). */
indices()136     const IndicesType& indices() const { return derived().indices(); }
137     /** \returns a reference to the stored array representing the permutation. */
indices()138     IndicesType& indices() { return derived().indices(); }
139 
140     /** Resizes to given size.
141       */
resize(Index newSize)142     inline void resize(Index newSize)
143     {
144       indices().resize(newSize);
145     }
146 
147     /** Sets *this to be the identity permutation matrix */
setIdentity()148     void setIdentity()
149     {
150       for(Index i = 0; i < size(); ++i)
151         indices().coeffRef(i) = i;
152     }
153 
154     /** Sets *this to be the identity permutation matrix of given size.
155       */
setIdentity(Index newSize)156     void setIdentity(Index newSize)
157     {
158       resize(newSize);
159       setIdentity();
160     }
161 
162     /** Multiplies *this by the transposition \f$(ij)\f$ on the left.
163       *
164       * \returns a reference to *this.
165       *
166       * \warning This is much slower than applyTranspositionOnTheRight(int,int):
167       * this has linear complexity and requires a lot of branching.
168       *
169       * \sa applyTranspositionOnTheRight(int,int)
170       */
applyTranspositionOnTheLeft(Index i,Index j)171     Derived& applyTranspositionOnTheLeft(Index i, Index j)
172     {
173       eigen_assert(i>=0 && j>=0 && i<size() && j<size());
174       for(Index k = 0; k < size(); ++k)
175       {
176         if(indices().coeff(k) == i) indices().coeffRef(k) = j;
177         else if(indices().coeff(k) == j) indices().coeffRef(k) = i;
178       }
179       return derived();
180     }
181 
182     /** Multiplies *this by the transposition \f$(ij)\f$ on the right.
183       *
184       * \returns a reference to *this.
185       *
186       * This is a fast operation, it only consists in swapping two indices.
187       *
188       * \sa applyTranspositionOnTheLeft(int,int)
189       */
applyTranspositionOnTheRight(Index i,Index j)190     Derived& applyTranspositionOnTheRight(Index i, Index j)
191     {
192       eigen_assert(i>=0 && j>=0 && i<size() && j<size());
193       std::swap(indices().coeffRef(i), indices().coeffRef(j));
194       return derived();
195     }
196 
197     /** \returns the inverse permutation matrix.
198       *
199       * \note \note_try_to_help_rvo
200       */
inverse()201     inline Transpose<PermutationBase> inverse() const
202     { return derived(); }
203     /** \returns the tranpose permutation matrix.
204       *
205       * \note \note_try_to_help_rvo
206       */
transpose()207     inline Transpose<PermutationBase> transpose() const
208     { return derived(); }
209 
210     /**** multiplication helpers to hopefully get RVO ****/
211 
212 
213 #ifndef EIGEN_PARSED_BY_DOXYGEN
214   protected:
215     template<typename OtherDerived>
assignTranspose(const PermutationBase<OtherDerived> & other)216     void assignTranspose(const PermutationBase<OtherDerived>& other)
217     {
218       for (int i=0; i<rows();++i) indices().coeffRef(other.indices().coeff(i)) = i;
219     }
220     template<typename Lhs,typename Rhs>
assignProduct(const Lhs & lhs,const Rhs & rhs)221     void assignProduct(const Lhs& lhs, const Rhs& rhs)
222     {
223       eigen_assert(lhs.cols() == rhs.rows());
224       for (int i=0; i<rows();++i) indices().coeffRef(i) = lhs.indices().coeff(rhs.indices().coeff(i));
225     }
226 #endif
227 
228   public:
229 
230     /** \returns the product permutation matrix.
231       *
232       * \note \note_try_to_help_rvo
233       */
234     template<typename Other>
235     inline PlainPermutationType operator*(const PermutationBase<Other>& other) const
236     { return PlainPermutationType(internal::PermPermProduct, derived(), other.derived()); }
237 
238     /** \returns the product of a permutation with another inverse permutation.
239       *
240       * \note \note_try_to_help_rvo
241       */
242     template<typename Other>
243     inline PlainPermutationType operator*(const Transpose<PermutationBase<Other> >& other) const
244     { return PlainPermutationType(internal::PermPermProduct, *this, other.eval()); }
245 
246     /** \returns the product of an inverse permutation with another permutation.
247       *
248       * \note \note_try_to_help_rvo
249       */
250     template<typename Other> friend
251     inline PlainPermutationType operator*(const Transpose<PermutationBase<Other> >& other, const PermutationBase& perm)
252     { return PlainPermutationType(internal::PermPermProduct, other.eval(), perm); }
253 
254     /** \returns the determinant of the permutation matrix, which is either 1 or -1 depending on the parity of the permutation.
255       *
256       * This function is O(\c n) procedure allocating a buffer of \c n booleans.
257       */
determinant()258     Index determinant() const
259     {
260       Index res = 1;
261       Index n = size();
262       Matrix<bool,RowsAtCompileTime,1,0,MaxRowsAtCompileTime> mask(n);
263       mask.fill(false);
264       Index r = 0;
265       while(r < n)
266       {
267         // search for the next seed
268         while(r<n && mask[r]) r++;
269         if(r>=n)
270           break;
271         // we got one, let's follow it until we are back to the seed
272         Index k0 = r++;
273         mask.coeffRef(k0) = true;
274         for(Index k=indices().coeff(k0); k!=k0; k=indices().coeff(k))
275         {
276           mask.coeffRef(k) = true;
277           res = -res;
278         }
279       }
280       return res;
281     }
282 
283   protected:
284 
285 };
286 
287 /** \class PermutationMatrix
288   * \ingroup Core_Module
289   *
290   * \brief Permutation matrix
291   *
292   * \param SizeAtCompileTime the number of rows/cols, or Dynamic
293   * \param MaxSizeAtCompileTime the maximum number of rows/cols, or Dynamic. This optional parameter defaults to SizeAtCompileTime. Most of the time, you should not have to specify it.
294   * \param IndexType the interger type of the indices
295   *
296   * This class represents a permutation matrix, internally stored as a vector of integers.
297   *
298   * \sa class PermutationBase, class PermutationWrapper, class DiagonalMatrix
299   */
300 
301 namespace internal {
302 template<int SizeAtCompileTime, int MaxSizeAtCompileTime, typename IndexType>
303 struct traits<PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime, IndexType> >
304  : traits<Matrix<IndexType,SizeAtCompileTime,SizeAtCompileTime,0,MaxSizeAtCompileTime,MaxSizeAtCompileTime> >
305 {
306   typedef IndexType Index;
307   typedef Matrix<IndexType, SizeAtCompileTime, 1, 0, MaxSizeAtCompileTime, 1> IndicesType;
308 };
309 }
310 
311 template<int SizeAtCompileTime, int MaxSizeAtCompileTime, typename IndexType>
312 class PermutationMatrix : public PermutationBase<PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime, IndexType> >
313 {
314     typedef PermutationBase<PermutationMatrix> Base;
315     typedef internal::traits<PermutationMatrix> Traits;
316   public:
317 
318     #ifndef EIGEN_PARSED_BY_DOXYGEN
319     typedef typename Traits::IndicesType IndicesType;
320     #endif
321 
322     inline PermutationMatrix()
323     {}
324 
325     /** Constructs an uninitialized permutation matrix of given size.
326       */
327     inline PermutationMatrix(int size) : m_indices(size)
328     {}
329 
330     /** Copy constructor. */
331     template<typename OtherDerived>
332     inline PermutationMatrix(const PermutationBase<OtherDerived>& other)
333       : m_indices(other.indices()) {}
334 
335     #ifndef EIGEN_PARSED_BY_DOXYGEN
336     /** Standard copy constructor. Defined only to prevent a default copy constructor
337       * from hiding the other templated constructor */
338     inline PermutationMatrix(const PermutationMatrix& other) : m_indices(other.indices()) {}
339     #endif
340 
341     /** Generic constructor from expression of the indices. The indices
342       * array has the meaning that the permutations sends each integer i to indices[i].
343       *
344       * \warning It is your responsibility to check that the indices array that you passes actually
345       * describes a permutation, i.e., each value between 0 and n-1 occurs exactly once, where n is the
346       * array's size.
347       */
348     template<typename Other>
349     explicit inline PermutationMatrix(const MatrixBase<Other>& a_indices) : m_indices(a_indices)
350     {}
351 
352     /** Convert the Transpositions \a tr to a permutation matrix */
353     template<typename Other>
354     explicit PermutationMatrix(const TranspositionsBase<Other>& tr)
355       : m_indices(tr.size())
356     {
357       *this = tr;
358     }
359 
360     /** Copies the other permutation into *this */
361     template<typename Other>
362     PermutationMatrix& operator=(const PermutationBase<Other>& other)
363     {
364       m_indices = other.indices();
365       return *this;
366     }
367 
368     /** Assignment from the Transpositions \a tr */
369     template<typename Other>
370     PermutationMatrix& operator=(const TranspositionsBase<Other>& tr)
371     {
372       return Base::operator=(tr.derived());
373     }
374 
375     #ifndef EIGEN_PARSED_BY_DOXYGEN
376     /** This is a special case of the templated operator=. Its purpose is to
377       * prevent a default operator= from hiding the templated operator=.
378       */
379     PermutationMatrix& operator=(const PermutationMatrix& other)
380     {
381       m_indices = other.m_indices;
382       return *this;
383     }
384     #endif
385 
386     /** const version of indices(). */
387     const IndicesType& indices() const { return m_indices; }
388     /** \returns a reference to the stored array representing the permutation. */
389     IndicesType& indices() { return m_indices; }
390 
391 
392     /**** multiplication helpers to hopefully get RVO ****/
393 
394 #ifndef EIGEN_PARSED_BY_DOXYGEN
395     template<typename Other>
396     PermutationMatrix(const Transpose<PermutationBase<Other> >& other)
397       : m_indices(other.nestedPermutation().size())
398     {
399       for (int i=0; i<m_indices.size();++i) m_indices.coeffRef(other.nestedPermutation().indices().coeff(i)) = i;
400     }
401     template<typename Lhs,typename Rhs>
402     PermutationMatrix(internal::PermPermProduct_t, const Lhs& lhs, const Rhs& rhs)
403       : m_indices(lhs.indices().size())
404     {
405       Base::assignProduct(lhs,rhs);
406     }
407 #endif
408 
409   protected:
410 
411     IndicesType m_indices;
412 };
413 
414 
415 namespace internal {
416 template<int SizeAtCompileTime, int MaxSizeAtCompileTime, typename IndexType, int _PacketAccess>
417 struct traits<Map<PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime, IndexType>,_PacketAccess> >
418  : traits<Matrix<IndexType,SizeAtCompileTime,SizeAtCompileTime,0,MaxSizeAtCompileTime,MaxSizeAtCompileTime> >
419 {
420   typedef IndexType Index;
421   typedef Map<const Matrix<IndexType, SizeAtCompileTime, 1, 0, MaxSizeAtCompileTime, 1>, _PacketAccess> IndicesType;
422 };
423 }
424 
425 template<int SizeAtCompileTime, int MaxSizeAtCompileTime, typename IndexType, int _PacketAccess>
426 class Map<PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime, IndexType>,_PacketAccess>
427   : public PermutationBase<Map<PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime, IndexType>,_PacketAccess> >
428 {
429     typedef PermutationBase<Map> Base;
430     typedef internal::traits<Map> Traits;
431   public:
432 
433     #ifndef EIGEN_PARSED_BY_DOXYGEN
434     typedef typename Traits::IndicesType IndicesType;
435     typedef typename IndicesType::Scalar Index;
436     #endif
437 
438     inline Map(const Index* indicesPtr)
439       : m_indices(indicesPtr)
440     {}
441 
442     inline Map(const Index* indicesPtr, Index size)
443       : m_indices(indicesPtr,size)
444     {}
445 
446     /** Copies the other permutation into *this */
447     template<typename Other>
448     Map& operator=(const PermutationBase<Other>& other)
449     { return Base::operator=(other.derived()); }
450 
451     /** Assignment from the Transpositions \a tr */
452     template<typename Other>
453     Map& operator=(const TranspositionsBase<Other>& tr)
454     { return Base::operator=(tr.derived()); }
455 
456     #ifndef EIGEN_PARSED_BY_DOXYGEN
457     /** This is a special case of the templated operator=. Its purpose is to
458       * prevent a default operator= from hiding the templated operator=.
459       */
460     Map& operator=(const Map& other)
461     {
462       m_indices = other.m_indices;
463       return *this;
464     }
465     #endif
466 
467     /** const version of indices(). */
468     const IndicesType& indices() const { return m_indices; }
469     /** \returns a reference to the stored array representing the permutation. */
470     IndicesType& indices() { return m_indices; }
471 
472   protected:
473 
474     IndicesType m_indices;
475 };
476 
477 /** \class PermutationWrapper
478   * \ingroup Core_Module
479   *
480   * \brief Class to view a vector of integers as a permutation matrix
481   *
482   * \param _IndicesType the type of the vector of integer (can be any compatible expression)
483   *
484   * This class allows to view any vector expression of integers as a permutation matrix.
485   *
486   * \sa class PermutationBase, class PermutationMatrix
487   */
488 
489 struct PermutationStorage {};
490 
491 template<typename _IndicesType> class TranspositionsWrapper;
492 namespace internal {
493 template<typename _IndicesType>
494 struct traits<PermutationWrapper<_IndicesType> >
495 {
496   typedef PermutationStorage StorageKind;
497   typedef typename _IndicesType::Scalar Scalar;
498   typedef typename _IndicesType::Scalar Index;
499   typedef _IndicesType IndicesType;
500   enum {
501     RowsAtCompileTime = _IndicesType::SizeAtCompileTime,
502     ColsAtCompileTime = _IndicesType::SizeAtCompileTime,
503     MaxRowsAtCompileTime = IndicesType::MaxRowsAtCompileTime,
504     MaxColsAtCompileTime = IndicesType::MaxColsAtCompileTime,
505     Flags = 0,
506     CoeffReadCost = _IndicesType::CoeffReadCost
507   };
508 };
509 }
510 
511 template<typename _IndicesType>
512 class PermutationWrapper : public PermutationBase<PermutationWrapper<_IndicesType> >
513 {
514     typedef PermutationBase<PermutationWrapper> Base;
515     typedef internal::traits<PermutationWrapper> Traits;
516   public:
517 
518     #ifndef EIGEN_PARSED_BY_DOXYGEN
519     typedef typename Traits::IndicesType IndicesType;
520     #endif
521 
522     inline PermutationWrapper(const IndicesType& a_indices)
523       : m_indices(a_indices)
524     {}
525 
526     /** const version of indices(). */
527     const typename internal::remove_all<typename IndicesType::Nested>::type&
528     indices() const { return m_indices; }
529 
530   protected:
531 
532     typename IndicesType::Nested m_indices;
533 };
534 
535 /** \returns the matrix with the permutation applied to the columns.
536   */
537 template<typename Derived, typename PermutationDerived>
538 inline const internal::permut_matrix_product_retval<PermutationDerived, Derived, OnTheRight>
539 operator*(const MatrixBase<Derived>& matrix,
540           const PermutationBase<PermutationDerived> &permutation)
541 {
542   return internal::permut_matrix_product_retval
543            <PermutationDerived, Derived, OnTheRight>
544            (permutation.derived(), matrix.derived());
545 }
546 
547 /** \returns the matrix with the permutation applied to the rows.
548   */
549 template<typename Derived, typename PermutationDerived>
550 inline const internal::permut_matrix_product_retval
551                <PermutationDerived, Derived, OnTheLeft>
552 operator*(const PermutationBase<PermutationDerived> &permutation,
553           const MatrixBase<Derived>& matrix)
554 {
555   return internal::permut_matrix_product_retval
556            <PermutationDerived, Derived, OnTheLeft>
557            (permutation.derived(), matrix.derived());
558 }
559 
560 namespace internal {
561 
562 template<typename PermutationType, typename MatrixType, int Side, bool Transposed>
563 struct traits<permut_matrix_product_retval<PermutationType, MatrixType, Side, Transposed> >
564 {
565   typedef typename MatrixType::PlainObject ReturnType;
566 };
567 
568 template<typename PermutationType, typename MatrixType, int Side, bool Transposed>
569 struct permut_matrix_product_retval
570  : public ReturnByValue<permut_matrix_product_retval<PermutationType, MatrixType, Side, Transposed> >
571 {
572     typedef typename remove_all<typename MatrixType::Nested>::type MatrixTypeNestedCleaned;
573     typedef typename MatrixType::Index Index;
574 
575     permut_matrix_product_retval(const PermutationType& perm, const MatrixType& matrix)
576       : m_permutation(perm), m_matrix(matrix)
577     {}
578 
579     inline Index rows() const { return m_matrix.rows(); }
580     inline Index cols() const { return m_matrix.cols(); }
581 
582     template<typename Dest> inline void evalTo(Dest& dst) const
583     {
584       const Index n = Side==OnTheLeft ? rows() : cols();
585       // FIXME we need an is_same for expression that is not sensitive to constness. For instance
586       // is_same_xpr<Block<const Matrix>, Block<Matrix> >::value should be true.
587       if(    is_same<MatrixTypeNestedCleaned,Dest>::value
588           && blas_traits<MatrixTypeNestedCleaned>::HasUsableDirectAccess
589           && blas_traits<Dest>::HasUsableDirectAccess
590           && extract_data(dst) == extract_data(m_matrix))
591       {
592         // apply the permutation inplace
593         Matrix<bool,PermutationType::RowsAtCompileTime,1,0,PermutationType::MaxRowsAtCompileTime> mask(m_permutation.size());
594         mask.fill(false);
595         Index r = 0;
596         while(r < m_permutation.size())
597         {
598           // search for the next seed
599           while(r<m_permutation.size() && mask[r]) r++;
600           if(r>=m_permutation.size())
601             break;
602           // we got one, let's follow it until we are back to the seed
603           Index k0 = r++;
604           Index kPrev = k0;
605           mask.coeffRef(k0) = true;
606           for(Index k=m_permutation.indices().coeff(k0); k!=k0; k=m_permutation.indices().coeff(k))
607           {
608                   Block<Dest, Side==OnTheLeft ? 1 : Dest::RowsAtCompileTime, Side==OnTheRight ? 1 : Dest::ColsAtCompileTime>(dst, k)
609             .swap(Block<Dest, Side==OnTheLeft ? 1 : Dest::RowsAtCompileTime, Side==OnTheRight ? 1 : Dest::ColsAtCompileTime>
610                        (dst,((Side==OnTheLeft) ^ Transposed) ? k0 : kPrev));
611 
612             mask.coeffRef(k) = true;
613             kPrev = k;
614           }
615         }
616       }
617       else
618       {
619         for(int i = 0; i < n; ++i)
620         {
621           Block<Dest, Side==OnTheLeft ? 1 : Dest::RowsAtCompileTime, Side==OnTheRight ? 1 : Dest::ColsAtCompileTime>
622                (dst, ((Side==OnTheLeft) ^ Transposed) ? m_permutation.indices().coeff(i) : i)
623 
624           =
625 
626           Block<const MatrixTypeNestedCleaned,Side==OnTheLeft ? 1 : MatrixType::RowsAtCompileTime,Side==OnTheRight ? 1 : MatrixType::ColsAtCompileTime>
627                (m_matrix, ((Side==OnTheRight) ^ Transposed) ? m_permutation.indices().coeff(i) : i);
628         }
629       }
630     }
631 
632   protected:
633     const PermutationType& m_permutation;
634     typename MatrixType::Nested m_matrix;
635 };
636 
637 /* Template partial specialization for transposed/inverse permutations */
638 
639 template<typename Derived>
640 struct traits<Transpose<PermutationBase<Derived> > >
641  : traits<Derived>
642 {};
643 
644 } // end namespace internal
645 
646 template<typename Derived>
647 class Transpose<PermutationBase<Derived> >
648   : public EigenBase<Transpose<PermutationBase<Derived> > >
649 {
650     typedef Derived PermutationType;
651     typedef typename PermutationType::IndicesType IndicesType;
652     typedef typename PermutationType::PlainPermutationType PlainPermutationType;
653   public:
654 
655     #ifndef EIGEN_PARSED_BY_DOXYGEN
656     typedef internal::traits<PermutationType> Traits;
657     typedef typename Derived::DenseMatrixType DenseMatrixType;
658     enum {
659       Flags = Traits::Flags,
660       CoeffReadCost = Traits::CoeffReadCost,
661       RowsAtCompileTime = Traits::RowsAtCompileTime,
662       ColsAtCompileTime = Traits::ColsAtCompileTime,
663       MaxRowsAtCompileTime = Traits::MaxRowsAtCompileTime,
664       MaxColsAtCompileTime = Traits::MaxColsAtCompileTime
665     };
666     typedef typename Traits::Scalar Scalar;
667     #endif
668 
669     Transpose(const PermutationType& p) : m_permutation(p) {}
670 
671     inline int rows() const { return m_permutation.rows(); }
672     inline int cols() const { return m_permutation.cols(); }
673 
674     #ifndef EIGEN_PARSED_BY_DOXYGEN
675     template<typename DenseDerived>
676     void evalTo(MatrixBase<DenseDerived>& other) const
677     {
678       other.setZero();
679       for (int i=0; i<rows();++i)
680         other.coeffRef(i, m_permutation.indices().coeff(i)) = typename DenseDerived::Scalar(1);
681     }
682     #endif
683 
684     /** \return the equivalent permutation matrix */
685     PlainPermutationType eval() const { return *this; }
686 
687     DenseMatrixType toDenseMatrix() const { return *this; }
688 
689     /** \returns the matrix with the inverse permutation applied to the columns.
690       */
691     template<typename OtherDerived> friend
692     inline const internal::permut_matrix_product_retval<PermutationType, OtherDerived, OnTheRight, true>
693     operator*(const MatrixBase<OtherDerived>& matrix, const Transpose& trPerm)
694     {
695       return internal::permut_matrix_product_retval<PermutationType, OtherDerived, OnTheRight, true>(trPerm.m_permutation, matrix.derived());
696     }
697 
698     /** \returns the matrix with the inverse permutation applied to the rows.
699       */
700     template<typename OtherDerived>
701     inline const internal::permut_matrix_product_retval<PermutationType, OtherDerived, OnTheLeft, true>
702     operator*(const MatrixBase<OtherDerived>& matrix) const
703     {
704       return internal::permut_matrix_product_retval<PermutationType, OtherDerived, OnTheLeft, true>(m_permutation, matrix.derived());
705     }
706 
707     const PermutationType& nestedPermutation() const { return m_permutation; }
708 
709   protected:
710     const PermutationType& m_permutation;
711 };
712 
713 template<typename Derived>
714 const PermutationWrapper<const Derived> MatrixBase<Derived>::asPermutation() const
715 {
716   return derived();
717 }
718 
719 } // end namespace Eigen
720 
721 #endif // EIGEN_PERMUTATIONMATRIX_H
722