mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-04-16 14:49:39 +08:00
* add PermutationMatrix
* DiagonalMatrix: - add MaxSizeAtCompileTime parameter - DiagonalOnTheLeft ---> OnTheLeft - fix bug in DiagonalMatrix::setIdentity()
This commit is contained in:
parent
1aaa9059c0
commit
955cd7f884
@ -158,6 +158,7 @@ ei_testing_print_summary()
|
||||
|
||||
message("")
|
||||
message("Configured Eigen ${EIGEN_VERSION_NUMBER}")
|
||||
message("")
|
||||
|
||||
string(TOLOWER "${CMAKE_GENERATOR}" cmake_generator_tolower)
|
||||
if(cmake_generator_tolower MATCHES "makefile")
|
||||
@ -166,14 +167,15 @@ if(cmake_generator_tolower MATCHES "makefile")
|
||||
message("Command | Description")
|
||||
message("--------------+----------------------------------------------------------------")
|
||||
message("make install | Install to ${CMAKE_INSTALL_PREFIX}")
|
||||
message(" | * To change that: cmake . -DCMAKE_INSTALL_PREFIX=yourpath")
|
||||
message(" | To change that: cmake . -DCMAKE_INSTALL_PREFIX=yourpath")
|
||||
message("make doc | Generate the API documentation, requires Doxygen & LaTeX")
|
||||
message("make check | Build and run the unit-tests")
|
||||
message("make check | Build and run the unit-tests. Read this page:")
|
||||
message(" | http://eigen.tuxfamily.org/index.php?title=Tests")
|
||||
message("make blas | Build BLAS library (not the same thing as Eigen)")
|
||||
message("--------------+----------------------------------------------------------------")
|
||||
else()
|
||||
message("To build/run the unit tests, read this page:")
|
||||
message(" http://eigen.tuxfamily.org/index.php?title=Tests")
|
||||
endif()
|
||||
|
||||
message("")
|
||||
message("To build/run the unit tests, read this page:")
|
||||
message(" http://eigen.tuxfamily.org/index.php?title=Tests")
|
||||
message("")
|
||||
|
@ -177,6 +177,7 @@ namespace Eigen {
|
||||
#include "src/Core/DiagonalMatrix.h"
|
||||
#include "src/Core/Diagonal.h"
|
||||
#include "src/Core/DiagonalProduct.h"
|
||||
#include "src/Core/PermutationMatrix.h"
|
||||
#include "src/Core/Redux.h"
|
||||
#include "src/Core/Visitor.h"
|
||||
#include "src/Core/Fuzzy.h"
|
||||
|
@ -66,7 +66,7 @@ class DiagonalBase : public AnyMatrixBase<Derived>
|
||||
inline int cols() const { return diagonal().size(); }
|
||||
|
||||
template<typename MatrixDerived>
|
||||
const DiagonalProduct<MatrixDerived, Derived, DiagonalOnTheLeft>
|
||||
const DiagonalProduct<MatrixDerived, Derived, OnTheLeft>
|
||||
operator*(const MatrixBase<MatrixDerived> &matrix) const;
|
||||
};
|
||||
|
||||
@ -88,16 +88,16 @@ void DiagonalBase<Derived>::evalTo(MatrixBase<DenseDerived> &other) const
|
||||
*
|
||||
* \sa class Matrix
|
||||
*/
|
||||
template<typename _Scalar, int _Size>
|
||||
struct ei_traits<DiagonalMatrix<_Scalar,_Size> >
|
||||
: ei_traits<Matrix<_Scalar,_Size,_Size> >
|
||||
template<typename _Scalar, int SizeAtCompileTime, int MaxSizeAtCompileTime>
|
||||
struct ei_traits<DiagonalMatrix<_Scalar,SizeAtCompileTime,MaxSizeAtCompileTime> >
|
||||
: ei_traits<Matrix<_Scalar,SizeAtCompileTime,SizeAtCompileTime,0,MaxSizeAtCompileTime,MaxSizeAtCompileTime> >
|
||||
{
|
||||
typedef Matrix<_Scalar,_Size,1> DiagonalVectorType;
|
||||
typedef Matrix<_Scalar,SizeAtCompileTime,1,0,MaxSizeAtCompileTime,1> DiagonalVectorType;
|
||||
};
|
||||
|
||||
template<typename _Scalar, int _Size>
|
||||
template<typename _Scalar, int SizeAtCompileTime, int MaxSizeAtCompileTime>
|
||||
class DiagonalMatrix
|
||||
: public DiagonalBase<DiagonalMatrix<_Scalar,_Size> >
|
||||
: public DiagonalBase<DiagonalMatrix<_Scalar,SizeAtCompileTime,MaxSizeAtCompileTime> >
|
||||
{
|
||||
public:
|
||||
|
||||
@ -156,8 +156,8 @@ class DiagonalMatrix
|
||||
inline void resize(int size) { m_diagonal.resize(size); }
|
||||
inline void setZero() { m_diagonal.setZero(); }
|
||||
inline void setZero(int size) { m_diagonal.setZero(size); }
|
||||
inline void setIdentity() { m_diagonal.setIdentity(); }
|
||||
inline void setIdentity(int size) { m_diagonal.setIdentity(size); }
|
||||
inline void setIdentity() { m_diagonal.setOnes(); }
|
||||
inline void setIdentity(int size) { m_diagonal.setOnes(size); }
|
||||
};
|
||||
|
||||
/** \class DiagonalWrapper
|
||||
|
@ -52,7 +52,7 @@ class DiagonalProduct : ei_no_assignment_operator,
|
||||
inline DiagonalProduct(const MatrixType& matrix, const DiagonalType& diagonal)
|
||||
: m_matrix(matrix), m_diagonal(diagonal)
|
||||
{
|
||||
ei_assert(diagonal.diagonal().size() == (ProductOrder == DiagonalOnTheLeft ? matrix.rows() : matrix.cols()));
|
||||
ei_assert(diagonal.diagonal().size() == (ProductOrder == OnTheLeft ? matrix.rows() : matrix.cols()));
|
||||
}
|
||||
|
||||
inline int rows() const { return m_matrix.rows(); }
|
||||
@ -60,7 +60,7 @@ class DiagonalProduct : ei_no_assignment_operator,
|
||||
|
||||
const Scalar coeff(int row, int col) const
|
||||
{
|
||||
return m_diagonal.diagonal().coeff(ProductOrder == DiagonalOnTheLeft ? row : col) * m_matrix.coeff(row, col);
|
||||
return m_diagonal.diagonal().coeff(ProductOrder == OnTheLeft ? row : col) * m_matrix.coeff(row, col);
|
||||
}
|
||||
|
||||
template<int LoadMode>
|
||||
@ -71,10 +71,10 @@ class DiagonalProduct : ei_no_assignment_operator,
|
||||
InnerSize = (MatrixType::Flags & RowMajorBit) ? MatrixType::ColsAtCompileTime : MatrixType::RowsAtCompileTime,
|
||||
DiagonalVectorPacketLoadMode = (LoadMode == Aligned && ((InnerSize%16) == 0)) ? Aligned : Unaligned
|
||||
};
|
||||
const int indexInDiagonalVector = ProductOrder == DiagonalOnTheLeft ? row : col;
|
||||
const int indexInDiagonalVector = ProductOrder == OnTheLeft ? row : col;
|
||||
|
||||
if((int(StorageOrder) == RowMajor && int(ProductOrder) == DiagonalOnTheLeft)
|
||||
||(int(StorageOrder) == ColMajor && int(ProductOrder) == DiagonalOnTheRight))
|
||||
if((int(StorageOrder) == RowMajor && int(ProductOrder) == OnTheLeft)
|
||||
||(int(StorageOrder) == ColMajor && int(ProductOrder) == OnTheRight))
|
||||
{
|
||||
return ei_pmul(m_matrix.template packet<LoadMode>(row, col),
|
||||
ei_pset1(m_diagonal.diagonal().coeff(indexInDiagonalVector)));
|
||||
@ -95,20 +95,20 @@ class DiagonalProduct : ei_no_assignment_operator,
|
||||
*/
|
||||
template<typename Derived>
|
||||
template<typename DiagonalDerived>
|
||||
inline const DiagonalProduct<Derived, DiagonalDerived, DiagonalOnTheRight>
|
||||
inline const DiagonalProduct<Derived, DiagonalDerived, OnTheRight>
|
||||
MatrixBase<Derived>::operator*(const DiagonalBase<DiagonalDerived> &diagonal) const
|
||||
{
|
||||
return DiagonalProduct<Derived, DiagonalDerived, DiagonalOnTheRight>(derived(), diagonal.derived());
|
||||
return DiagonalProduct<Derived, DiagonalDerived, OnTheRight>(derived(), diagonal.derived());
|
||||
}
|
||||
|
||||
/** \returns the diagonal matrix product of \c *this by the matrix \a matrix.
|
||||
*/
|
||||
template<typename DiagonalDerived>
|
||||
template<typename MatrixDerived>
|
||||
inline const DiagonalProduct<MatrixDerived, DiagonalDerived, DiagonalOnTheLeft>
|
||||
inline const DiagonalProduct<MatrixDerived, DiagonalDerived, OnTheLeft>
|
||||
DiagonalBase<DiagonalDerived>::operator*(const MatrixBase<MatrixDerived> &matrix) const
|
||||
{
|
||||
return DiagonalProduct<MatrixDerived, DiagonalDerived, DiagonalOnTheLeft>(matrix.derived(), derived());
|
||||
return DiagonalProduct<MatrixDerived, DiagonalDerived, OnTheLeft>(matrix.derived(), derived());
|
||||
}
|
||||
|
||||
|
||||
|
@ -410,7 +410,7 @@ template<typename Derived> class MatrixBase
|
||||
void applyOnTheRight(const AnyMatrixBase<OtherDerived>& other);
|
||||
|
||||
template<typename DiagonalDerived>
|
||||
const DiagonalProduct<Derived, DiagonalDerived, DiagonalOnTheRight>
|
||||
const DiagonalProduct<Derived, DiagonalDerived, OnTheRight>
|
||||
operator*(const DiagonalBase<DiagonalDerived> &diagonal) const;
|
||||
|
||||
template<typename OtherDerived>
|
||||
|
205
Eigen/src/Core/PermutationMatrix.h
Normal file
205
Eigen/src/Core/PermutationMatrix.h
Normal file
@ -0,0 +1,205 @@
|
||||
// This file is part of Eigen, a lightweight C++ template library
|
||||
// for linear algebra.
|
||||
//
|
||||
// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
|
||||
//
|
||||
// Eigen is free software; you can redistribute it and/or
|
||||
// modify it under the terms of the GNU Lesser General Public
|
||||
// License as published by the Free Software Foundation; either
|
||||
// version 3 of the License, or (at your option) any later version.
|
||||
//
|
||||
// Alternatively, you can redistribute it and/or
|
||||
// modify it under the terms of the GNU General Public License as
|
||||
// published by the Free Software Foundation; either version 2 of
|
||||
// the License, or (at your option) any later version.
|
||||
//
|
||||
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
|
||||
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
|
||||
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
|
||||
// GNU General Public License for more details.
|
||||
//
|
||||
// You should have received a copy of the GNU Lesser General Public
|
||||
// License and a copy of the GNU General Public License along with
|
||||
// Eigen. If not, see <http://www.gnu.org/licenses/>.
|
||||
|
||||
#ifndef EIGEN_PERMUTATIONMATRIX_H
|
||||
#define EIGEN_PERMUTATIONMATRIX_H
|
||||
|
||||
/** \nonstableyet
|
||||
* \class PermutationMatrix
|
||||
*
|
||||
* \brief Permutation matrix
|
||||
*
|
||||
* \param SizeAtCompileTime the number of rows/cols, or Dynamic
|
||||
* \param MaxSizeAtCompileTime the maximum number of rows/cols, or Dynamic. This optional parameter defaults to SizeAtCompileTime.
|
||||
*
|
||||
* This class represents a permutation matrix, internally stored as a vector of integers.
|
||||
* The convention followed here is the same as on <a href="http://en.wikipedia.org/wiki/Permutation_matrix">Wikipedia</a>,
|
||||
* namely: the matrix of permutation \a p is the matrix such that on each row \a i, the only nonzero coefficient is
|
||||
* in column p(i).
|
||||
*
|
||||
* \sa class DiagonalMatrix
|
||||
*/
|
||||
template<int SizeAtCompileTime, int MaxSizeAtCompileTime = SizeAtCompileTime> class PermutationMatrix;
|
||||
template<typename PermutationType, typename MatrixType, int Side> struct ei_permut_matrix_product_retval;
|
||||
|
||||
template<int SizeAtCompileTime, int MaxSizeAtCompileTime>
|
||||
struct ei_traits<PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime> >
|
||||
: ei_traits<Matrix<int,SizeAtCompileTime,SizeAtCompileTime,0,MaxSizeAtCompileTime,MaxSizeAtCompileTime> >
|
||||
{};
|
||||
|
||||
template<int SizeAtCompileTime, int MaxSizeAtCompileTime>
|
||||
class PermutationMatrix : public AnyMatrixBase<PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime> >
|
||||
{
|
||||
public:
|
||||
|
||||
typedef ei_traits<PermutationMatrix> Traits;
|
||||
typedef Matrix<int,SizeAtCompileTime,SizeAtCompileTime,0,MaxSizeAtCompileTime,MaxSizeAtCompileTime>
|
||||
DenseMatrixType;
|
||||
enum {
|
||||
Flags = Traits::Flags,
|
||||
CoeffReadCost = Traits::CoeffReadCost,
|
||||
RowsAtCompileTime = Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = Traits::ColsAtCompileTime,
|
||||
MaxRowsAtCompileTime = Traits::MaxRowsAtCompileTime,
|
||||
MaxColsAtCompileTime = Traits::MaxColsAtCompileTime
|
||||
};
|
||||
typedef typename Traits::Scalar Scalar;
|
||||
|
||||
typedef Matrix<int, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1> IndicesType;
|
||||
|
||||
inline PermutationMatrix()
|
||||
{
|
||||
}
|
||||
|
||||
template<int OtherSize, int OtherMaxSize>
|
||||
inline PermutationMatrix(const PermutationMatrix<OtherSize, OtherMaxSize>& other)
|
||||
: m_indices(other.indices()) {}
|
||||
|
||||
/** copy constructor. prevent a default copy constructor from hiding the other templated constructor */
|
||||
inline PermutationMatrix(const PermutationMatrix& other) : m_indices(other.indices()) {}
|
||||
|
||||
/** generic constructor from expression of the indices */
|
||||
template<typename Other>
|
||||
explicit inline PermutationMatrix(const MatrixBase<Other>& other) : m_indices(other)
|
||||
{}
|
||||
|
||||
template<int OtherSize, int OtherMaxSize>
|
||||
PermutationMatrix& operator=(const PermutationMatrix<OtherSize, OtherMaxSize>& other)
|
||||
{
|
||||
m_indices = other.indices();
|
||||
return *this;
|
||||
}
|
||||
|
||||
/** This is a special case of the templated operator=. Its purpose is to
|
||||
* prevent a default operator= from hiding the templated operator=.
|
||||
*/
|
||||
PermutationMatrix& operator=(const PermutationMatrix& other)
|
||||
{
|
||||
m_indices = other.m_indices();
|
||||
return *this;
|
||||
}
|
||||
|
||||
inline PermutationMatrix(int rows, int cols) : m_indices(rows)
|
||||
{
|
||||
ei_assert(rows == cols);
|
||||
}
|
||||
|
||||
/** \returns the number of columns */
|
||||
inline int rows() const { return m_indices.size(); }
|
||||
|
||||
/** \returns the number of rows */
|
||||
inline int cols() const { return m_indices.size(); }
|
||||
|
||||
template<typename DenseDerived>
|
||||
void evalTo(MatrixBase<DenseDerived>& other) const
|
||||
{
|
||||
other.setZero();
|
||||
for (int i=0; i<rows();++i)
|
||||
other.coeffRef(i,m_indices.coeff(i)) = typename DenseDerived::Scalar(1);
|
||||
}
|
||||
|
||||
DenseMatrixType toDenseMatrix() const
|
||||
{
|
||||
return *this;
|
||||
}
|
||||
|
||||
const IndicesType& indices() const { return m_indices; }
|
||||
IndicesType& indices() { return m_indices; }
|
||||
|
||||
protected:
|
||||
|
||||
IndicesType m_indices;
|
||||
};
|
||||
|
||||
/** \returns the matrix with the permutation applied to the columns.
|
||||
*/
|
||||
template<typename Derived, int SizeAtCompileTime, int MaxSizeAtCompileTime>
|
||||
inline const ei_permut_matrix_product_retval<PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime>, Derived, OnTheRight>
|
||||
operator*(const MatrixBase<Derived>& matrix,
|
||||
const PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime> &permutation)
|
||||
{
|
||||
return ei_permut_matrix_product_retval
|
||||
<PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime>, Derived, OnTheRight>
|
||||
(permutation, matrix.derived());
|
||||
}
|
||||
|
||||
/** \returns the matrix with the permutation applied to the rows.
|
||||
*/
|
||||
template<typename Derived, int SizeAtCompileTime, int MaxSizeAtCompileTime>
|
||||
inline const ei_permut_matrix_product_retval
|
||||
<PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime>, Derived, OnTheLeft>
|
||||
operator*(const PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime> &permutation,
|
||||
const MatrixBase<Derived>& matrix)
|
||||
{
|
||||
return ei_permut_matrix_product_retval
|
||||
<PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime>, Derived, OnTheLeft>
|
||||
(permutation, matrix.derived());
|
||||
}
|
||||
|
||||
template<typename PermutationType, typename MatrixType, int Side>
|
||||
struct ei_traits<ei_permut_matrix_product_retval<PermutationType, MatrixType, Side> >
|
||||
{
|
||||
typedef typename MatrixType::PlainMatrixType ReturnMatrixType;
|
||||
};
|
||||
|
||||
template<typename PermutationType, typename MatrixType, int Side>
|
||||
struct ei_permut_matrix_product_retval
|
||||
: public ReturnByValue<ei_permut_matrix_product_retval<PermutationType, MatrixType, Side> >
|
||||
{
|
||||
typedef typename ei_cleantype<typename MatrixType::Nested>::type MatrixTypeNestedCleaned;
|
||||
|
||||
ei_permut_matrix_product_retval(const PermutationType& perm, const MatrixType& matrix)
|
||||
: m_permutation(perm), m_matrix(matrix)
|
||||
{}
|
||||
|
||||
inline int rows() const { return m_matrix.rows(); }
|
||||
inline int cols() const { return m_matrix.cols(); }
|
||||
|
||||
template<typename Dest> inline void evalTo(Dest& dst) const
|
||||
{
|
||||
const int n = Side==OnTheLeft ? rows() : cols();
|
||||
for(int i = 0; i < n; ++i)
|
||||
{
|
||||
Block<
|
||||
Dest,
|
||||
Side==OnTheLeft ? 1 : Dest::RowsAtCompileTime,
|
||||
Side==OnTheRight ? 1 : Dest::ColsAtCompileTime
|
||||
>(dst, Side==OnTheRight ? m_permutation.indices().coeff(i) : i)
|
||||
|
||||
=
|
||||
|
||||
Block<
|
||||
MatrixTypeNestedCleaned,
|
||||
Side==OnTheLeft ? 1 : MatrixType::RowsAtCompileTime,
|
||||
Side==OnTheRight ? 1 : MatrixType::ColsAtCompileTime
|
||||
>(m_matrix, Side==OnTheLeft ? m_permutation.indices().coeff(i) : i);
|
||||
}
|
||||
}
|
||||
|
||||
protected:
|
||||
const PermutationType& m_permutation;
|
||||
const typename MatrixType::Nested m_matrix;
|
||||
};
|
||||
|
||||
#endif // EIGEN_PERMUTATIONMATRIX_H
|
@ -194,8 +194,6 @@ const unsigned int SelfAdjoint = SelfAdjointBit;
|
||||
const unsigned int UnitUpperTriangular = UpperTriangularBit | UnitDiagBit;
|
||||
const unsigned int UnitLowerTriangular = LowerTriangularBit | UnitDiagBit;
|
||||
|
||||
enum { DiagonalOnTheLeft, DiagonalOnTheRight };
|
||||
|
||||
enum { Unaligned=0, Aligned=1 };
|
||||
enum { AsRequested=0, EnforceAlignedAccess=2 };
|
||||
enum { ConditionalJumpCost = 5 };
|
||||
@ -232,7 +230,6 @@ enum {
|
||||
DontAlign = 0x2
|
||||
};
|
||||
|
||||
// used for the solvers
|
||||
enum {
|
||||
OnTheLeft = 1,
|
||||
OnTheRight = 2
|
||||
|
@ -53,7 +53,7 @@ template<typename Derived, typename Lhs, typename Rhs> class ProductBase;
|
||||
|
||||
template<typename Derived> class DiagonalBase;
|
||||
template<typename _DiagonalVectorType> class DiagonalWrapper;
|
||||
template<typename _Scalar, int _Size> class DiagonalMatrix;
|
||||
template<typename _Scalar, int SizeAtCompileTime, int MaxSizeAtCompileTime=SizeAtCompileTime> class DiagonalMatrix;
|
||||
template<typename MatrixType, typename DiagonalType, int ProductOrder> class DiagonalProduct;
|
||||
template<typename MatrixType, int Index> class Diagonal;
|
||||
|
||||
|
@ -3,7 +3,7 @@
|
||||
|
||||
if [ $# == 0 -o $# -ge 3 ]
|
||||
then
|
||||
echo "usage: ./mctestr regexp [jobs]"
|
||||
echo "usage: ./check regexp [jobs]"
|
||||
echo " makes and runs tests matching the regexp, with <jobs> concurrent make jobs"
|
||||
exit 0
|
||||
fi
|
||||
|
@ -158,6 +158,7 @@ ei_add_test(umeyama)
|
||||
ei_add_test(householder)
|
||||
ei_add_test(swap)
|
||||
ei_add_test(conservative_resize)
|
||||
ei_add_test(permutationmatrices)
|
||||
|
||||
ei_add_property(EIGEN_TESTING_SUMMARY "CXX: ${CMAKE_CXX_COMPILER}\n")
|
||||
if(CMAKE_COMPILER_IS_GNUCXX)
|
||||
|
@ -2,7 +2,7 @@
|
||||
|
||||
if [ $# == 0 -o $# -ge 3 ]
|
||||
then
|
||||
echo "usage: ./maketests regexp [jobs]"
|
||||
echo "usage: ./buildtests regexp [jobs]"
|
||||
echo " makes tests matching the regexp, with <jobs> concurrent make jobs"
|
||||
exit 0
|
||||
fi
|
||||
|
89
test/permutationmatrices.cpp
Normal file
89
test/permutationmatrices.cpp
Normal file
@ -0,0 +1,89 @@
|
||||
// This file is part of Eigen, a lightweight C++ template library
|
||||
// for linear algebra.
|
||||
//
|
||||
// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
|
||||
//
|
||||
// Eigen is free software; you can redistribute it and/or
|
||||
// modify it under the terms of the GNU Lesser General Public
|
||||
// License as published by the Free Software Foundation; either
|
||||
// version 3 of the License, or (at your option) any later version.
|
||||
//
|
||||
// Alternatively, you can redistribute it and/or
|
||||
// modify it under the terms of the GNU General Public License as
|
||||
// published by the Free Software Foundation; either version 2 of
|
||||
// the License, or (at your option) any later version.
|
||||
//
|
||||
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
|
||||
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
|
||||
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
|
||||
// GNU General Public License for more details.
|
||||
//
|
||||
// You should have received a copy of the GNU Lesser General Public
|
||||
// License and a copy of the GNU General Public License along with
|
||||
// Eigen. If not, see <http://www.gnu.org/licenses/>.
|
||||
|
||||
#include "main.h"
|
||||
|
||||
template<typename PermutationVectorType>
|
||||
void randomPermutationVector(PermutationVectorType& v, int size)
|
||||
{
|
||||
typedef typename PermutationVectorType::Scalar Scalar;
|
||||
v.resize(size);
|
||||
for(int i = 0; i < size; ++i) v(i) = Scalar(i);
|
||||
if(size == 1) return;
|
||||
for(int n = 0; n < 3 * size; ++n)
|
||||
{
|
||||
int i = ei_random<int>(0, size-1);
|
||||
int j;
|
||||
do j = ei_random<int>(0, size-1); while(j==i);
|
||||
std::swap(v(i), v(j));
|
||||
}
|
||||
}
|
||||
|
||||
using namespace std;
|
||||
template<typename MatrixType> void permutationmatrices(const MatrixType& m)
|
||||
{
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef typename MatrixType::RealScalar RealScalar;
|
||||
enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime,
|
||||
Options = MatrixType::Options };
|
||||
typedef PermutationMatrix<Rows> LeftPermutationType;
|
||||
typedef Matrix<int, Rows, 1> LeftPermutationVectorType;
|
||||
typedef PermutationMatrix<Cols> RightPermutationType;
|
||||
typedef Matrix<int, Cols, 1> RightPermutationVectorType;
|
||||
|
||||
int rows = m.rows();
|
||||
int cols = m.cols();
|
||||
|
||||
MatrixType m_original = MatrixType::Random(rows,cols);
|
||||
LeftPermutationVectorType lv;
|
||||
randomPermutationVector(lv, rows);
|
||||
LeftPermutationType lp(lv);
|
||||
RightPermutationVectorType rv;
|
||||
randomPermutationVector(rv, cols);
|
||||
RightPermutationType rp(rv);
|
||||
MatrixType m_permuted = lp * m_original * rp;
|
||||
|
||||
int i = ei_random<int>(0, rows-1);
|
||||
int j = ei_random<int>(0, cols-1);
|
||||
|
||||
VERIFY_IS_APPROX(m_original(lv(i),j), m_permuted(i,rv(j)));
|
||||
|
||||
Matrix<Scalar,Rows,Rows> lm(lp);
|
||||
Matrix<Scalar,Cols,Cols> rm(rp);
|
||||
|
||||
VERIFY_IS_APPROX(m_permuted, lm*m_original*rm);
|
||||
}
|
||||
|
||||
void test_permutationmatrices()
|
||||
{
|
||||
for(int i = 0; i < g_repeat; i++) {
|
||||
CALL_SUBTEST_1( permutationmatrices(Matrix<float, 1, 1>()) );
|
||||
CALL_SUBTEST_2( permutationmatrices(Matrix3f()) );
|
||||
CALL_SUBTEST_3( permutationmatrices(Matrix<double,3,3,RowMajor>()) );
|
||||
CALL_SUBTEST_4( permutationmatrices(Matrix4d()) );
|
||||
CALL_SUBTEST_5( permutationmatrices(Matrix<double,4,6>()) );
|
||||
CALL_SUBTEST_6( permutationmatrices(Matrix<double,Dynamic,Dynamic,RowMajor>(20, 30)) );
|
||||
CALL_SUBTEST_7( permutationmatrices(MatrixXcf(15, 10)) );
|
||||
}
|
||||
}
|
Loading…
x
Reference in New Issue
Block a user