PermutationWrapper< _IndicesType > Class Template Reference

Class to view a vector of integers as a permutation matrix. More...

#include <PermutationMatrix.h>

+ Inheritance diagram for PermutationWrapper< _IndicesType >:

Public Member Functions

PermutationWrapper< _IndicesType > & applyTranspositionOnTheLeft (Index i, Index j)
 
PermutationWrapper< _IndicesType > & applyTranspositionOnTheRight (Index i, Index j)
 
Index cols () const
 
PermutationWrapper< _IndicesType > & derived ()
 
const PermutationWrapper< _IndicesType > & derived () const
 
IndicesType & indices ()
 
const internal::remove_all< typenameIndicesType::Nested >::type & indices () const
 
Transpose< PermutationBaseinverse () const
 
PlainPermutationType operator* (const PermutationBase< Other > &other) const
 
PlainPermutationType operator* (const Transpose< PermutationBase< Other > > &other) const
 
void resize (Index size)
 
Index rows () const
 
void setIdentity ()
 
void setIdentity (Index size)
 
Index size () const
 
DenseMatrixType toDenseMatrix () const
 
Transpose< PermutationBasetranspose () const
 

Detailed Description

template<typename _IndicesType>
class Eigen::PermutationWrapper< _IndicesType >

Class to view a vector of integers as a permutation matrix.

Parameters
_IndicesTypethe type of the vector of integer (can be any compatible expression)

This class allows to view any vector expression of integers as a permutation matrix.

See also
class PermutationBase, class PermutationMatrix

Member Function Documentation

◆ applyTranspositionOnTheLeft()

PermutationWrapper< _IndicesType > & applyTranspositionOnTheLeft ( Index i,
Index j )
inlineinherited

Multiplies *this by the transposition $(ij)$ on the left.

Returns
a reference to *this.
Warning
This is much slower than applyTranspositionOnTheRight(int,int): this has linear complexity and requires a lot of branching.
See also
applyTranspositionOnTheRight(int,int)

◆ applyTranspositionOnTheRight()

PermutationWrapper< _IndicesType > & applyTranspositionOnTheRight ( Index i,
Index j )
inlineinherited

Multiplies *this by the transposition $(ij)$ on the right.

Returns
a reference to *this.

This is a fast operation, it only consists in swapping two indices.

See also
applyTranspositionOnTheLeft(int,int)

◆ cols()

Index cols ( ) const
inlineinherited
Returns
the number of columns

◆ derived() [1/2]

PermutationWrapper< _IndicesType > & derived ( )
inlineinherited
Returns
a reference to the derived object

◆ derived() [2/2]

const PermutationWrapper< _IndicesType > & derived ( ) const
inlineinherited
Returns
a const reference to the derived object

◆ indices() [1/2]

IndicesType & indices ( )
inlineinherited
Returns
a reference to the stored array representing the permutation.

◆ indices() [2/2]

template<typename _IndicesType>
const internal::remove_all< typenameIndicesType::Nested >::type & indices ( ) const
inline

const version of indices().

◆ inverse()

Transpose< PermutationBase > inverse ( ) const
inlineinherited
Returns
the inverse permutation matrix.
Note
This function returns the result by value. In order to make that efficient, it is implemented as just a return statement using a special constructor, hopefully allowing the compiler to perform a RVO (return value optimization).

◆ operator*() [1/2]

PlainPermutationType operator* ( const PermutationBase< Other > & other) const
inlineinherited
Returns
the product permutation matrix.
Note
This function returns the result by value. In order to make that efficient, it is implemented as just a return statement using a special constructor, hopefully allowing the compiler to perform a RVO (return value optimization).

◆ operator*() [2/2]

PlainPermutationType operator* ( const Transpose< PermutationBase< Other > > & other) const
inlineinherited
Returns
the product of a permutation with another inverse permutation.
Note
This function returns the result by value. In order to make that efficient, it is implemented as just a return statement using a special constructor, hopefully allowing the compiler to perform a RVO (return value optimization).

◆ resize()

void resize ( Index size)
inlineinherited

Resizes to given size.

◆ rows()

Index rows ( ) const
inlineinherited
Returns
the number of rows

◆ setIdentity() [1/2]

void setIdentity ( )
inlineinherited

Sets *this to be the identity permutation matrix

◆ setIdentity() [2/2]

void setIdentity ( Index size)
inlineinherited

Sets *this to be the identity permutation matrix of given size.

◆ size()

Index size ( ) const
inlineinherited
Returns
the size of a side of the respective square matrix, i.e., the number of indices

◆ toDenseMatrix()

DenseMatrixType toDenseMatrix ( ) const
inlineinherited
Returns
a Matrix object initialized from this permutation matrix. Notice that it is inefficient to return this Matrix object by value. For efficiency, favor using the Matrix constructor taking EigenBase objects.

◆ transpose()

Transpose< PermutationBase > transpose ( ) const
inlineinherited
Returns
the tranpose permutation matrix.
Note
This function returns the result by value. In order to make that efficient, it is implemented as just a return statement using a special constructor, hopefully allowing the compiler to perform a RVO (return value optimization).

The documentation for this class was generated from the following file: