Propagate all five matrix template parameters to members and temporaries of decomposition classes. One particular advantage of this is that decomposing matrices with max sizes known at compile time will not allocate.

NOTE: The ComplexEigenSolver class currently _does_ allocate (line 135 of Eigenvalues/ComplexEigenSolver.h), but the reason appears to be in the implementation of matrix-matrix products, and not in the decomposition itself.
The nomalloc unit test has been extended to verify that decompositions do not allocate when max sizes are specified. There are currently two workarounds to prevent the test from failing (see comments in test/nomalloc.cpp), both of which are related to matrix products that allocate on the stack.
This commit is contained in:
Adolfo Rodriguez Tsouroukdissian 2010-03-08 19:31:27 +01:00
parent 898762529e
commit 5a36f4a8d1
15 changed files with 192 additions and 69 deletions

View File

@ -56,11 +56,18 @@ template<typename _MatrixType> class LDLT
{
public:
typedef _MatrixType MatrixType;
enum {
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
Options = MatrixType::Options,
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
};
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> IntColVectorType;
typedef Matrix<int, 1, MatrixType::RowsAtCompileTime> IntRowVectorType;
typedef Matrix<Scalar, ColsAtCompileTime, 1, Options, MaxColsAtCompileTime, 1> VectorType;
typedef Matrix<int, RowsAtCompileTime, 1, Options, MaxRowsAtCompileTime, 1> IntColVectorType;
typedef Matrix<int, 1, RowsAtCompileTime, Options, 1, MaxRowsAtCompileTime> IntRowVectorType;
/** \brief Default Constructor.
*
@ -201,7 +208,7 @@ LDLT<MatrixType>& LDLT<MatrixType>::compute(const MatrixType& a)
// By using a temorary, packet-aligned products are guarenteed. In the LLT
// case this is unnecessary because the diagonal is included and will always
// have optimal alignment.
Matrix<Scalar,MatrixType::RowsAtCompileTime,1> _temporary(size);
Matrix<Scalar, RowsAtCompileTime, 1, Options, MaxRowsAtCompileTime, 1> _temporary(size);
for (int j = 0; j < size; ++j)
{

View File

@ -57,9 +57,15 @@ template<typename _MatrixType, int _UpLo> class LLT
{
public:
typedef _MatrixType MatrixType;
enum {
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
Options = MatrixType::Options,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
};
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
typedef Matrix<Scalar, ColsAtCompileTime, 1, Options, MaxColsAtCompileTime, 1> VectorType;
enum {
PacketSize = ei_packet_traits<Scalar>::size,

View File

@ -41,11 +41,18 @@ template<typename _MatrixType> class ComplexEigenSolver
{
public:
typedef _MatrixType MatrixType;
enum {
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
Options = MatrixType::Options,
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
};
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef std::complex<RealScalar> Complex;
typedef Matrix<Complex, MatrixType::ColsAtCompileTime,1> EigenvalueType;
typedef Matrix<Complex, MatrixType::RowsAtCompileTime,MatrixType::ColsAtCompileTime> EigenvectorType;
typedef Matrix<Complex, ColsAtCompileTime, 1, Options, MaxColsAtCompileTime, 1> EigenvalueType;
typedef Matrix<Complex, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, ColsAtCompileTime> EigenvectorType;
/**
* \brief Default Constructor.

View File

@ -45,10 +45,17 @@ template<typename _MatrixType> class ComplexSchur
{
public:
typedef _MatrixType MatrixType;
enum {
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
Options = MatrixType::Options,
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
};
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef std::complex<RealScalar> Complex;
typedef Matrix<Complex, MatrixType::RowsAtCompileTime,MatrixType::ColsAtCompileTime> ComplexMatrixType;
typedef Matrix<Complex, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> ComplexMatrixType;
enum {
Size = MatrixType::RowsAtCompileTime
};

View File

@ -45,13 +45,19 @@ template<typename _MatrixType> class EigenSolver
public:
typedef _MatrixType MatrixType;
enum {
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
Options = MatrixType::Options,
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
};
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef std::complex<RealScalar> Complex;
typedef Matrix<Complex, MatrixType::ColsAtCompileTime, 1> EigenvalueType;
typedef Matrix<Complex, MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime> EigenvectorType;
typedef Matrix<RealScalar, MatrixType::ColsAtCompileTime, 1> RealVectorType;
typedef Matrix<RealScalar, Dynamic, 1> RealVectorTypeX;
typedef Matrix<Complex, ColsAtCompileTime, 1, Options, MaxColsAtCompileTime, 1> EigenvalueType;
typedef Matrix<Complex, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> EigenvectorType;
typedef Matrix<RealScalar, ColsAtCompileTime, 1, Options, MaxColsAtCompileTime, 1> RealVectorType;
/**
* \brief Default Constructor.

View File

@ -44,17 +44,17 @@ template<typename _MatrixType> class HessenbergDecomposition
public:
typedef _MatrixType MatrixType;
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
enum {
Size = MatrixType::RowsAtCompileTime,
SizeMinusOne = MatrixType::RowsAtCompileTime==Dynamic
? Dynamic
: MatrixType::RowsAtCompileTime-1
SizeMinusOne = Size == Dynamic ? Dynamic : Size - 1,
Options = MatrixType::Options,
MaxSize = MatrixType::MaxRowsAtCompileTime,
MaxSizeMinusOne = MaxSize == Dynamic ? Dynamic : MaxSize - 1
};
typedef Matrix<Scalar, SizeMinusOne, 1> CoeffVectorType;
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef Matrix<Scalar, SizeMinusOne, 1, Options, MaxSizeMinusOne, 1> CoeffVectorType;
typedef Matrix<Scalar, 1, Size, Options, 1, MaxSize> VectorType;
/** This constructor initializes a HessenbergDecomposition object for
* further use with HessenbergDecomposition::compute()
@ -143,7 +143,7 @@ void HessenbergDecomposition<MatrixType>::_compute(MatrixType& matA, CoeffVector
{
assert(matA.rows()==matA.cols());
int n = matA.rows();
Matrix<Scalar,1,Dynamic> temp(n);
VectorType temp(n);
for (int i = 0; i<n-1; ++i)
{
// let's consider the vector v = i-th column starting at position i+1
@ -174,7 +174,7 @@ HessenbergDecomposition<MatrixType>::matrixQ() const
{
int n = m_matrix.rows();
MatrixType matQ = MatrixType::Identity(n,n);
Matrix<Scalar,1,MatrixType::ColsAtCompileTime> temp(n);
VectorType temp(n);
for (int i = n-2; i>=0; i--)
{
matQ.corner(BottomRight,n-i-1,n-i-1)

View File

@ -42,13 +42,17 @@ template<typename _MatrixType> class SelfAdjointEigenSolver
{
public:
enum {Size = _MatrixType::RowsAtCompileTime };
typedef _MatrixType MatrixType;
enum {
Size = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
Options = MatrixType::Options,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
};
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef std::complex<RealScalar> Complex;
typedef Matrix<RealScalar, MatrixType::ColsAtCompileTime, 1> RealVectorType;
typedef Matrix<RealScalar, Dynamic, 1> RealVectorTypeX;
typedef Matrix<RealScalar, ColsAtCompileTime, 1, Options, MaxColsAtCompileTime, 1> RealVectorType;
typedef Tridiagonalization<MatrixType> TridiagonalizationType;
// typedef typename TridiagonalizationType::TridiagonalMatrixType TridiagonalMatrixType;

View File

@ -50,16 +50,18 @@ template<typename _MatrixType> class Tridiagonalization
enum {
Size = MatrixType::RowsAtCompileTime,
SizeMinusOne = MatrixType::RowsAtCompileTime==Dynamic
? Dynamic
: MatrixType::RowsAtCompileTime-1,
SizeMinusOne = Size == Dynamic ? Dynamic : Size - 1,
Options = MatrixType::Options,
MaxSize = MatrixType::MaxRowsAtCompileTime,
MaxSizeMinusOne = MaxSize == Dynamic ? Dynamic : MaxSize - 1,
PacketSize = ei_packet_traits<Scalar>::size
};
typedef Matrix<Scalar, SizeMinusOne, 1> CoeffVectorType;
typedef Matrix<RealScalar, Size, 1> DiagonalType;
typedef Matrix<RealScalar, SizeMinusOne, 1> SubDiagonalType;
typedef Matrix<Scalar, SizeMinusOne, 1, Options, MaxSizeMinusOne, 1> CoeffVectorType;
typedef Matrix<RealScalar, Size, 1, Options, MaxSize, 1> DiagonalType;
typedef Matrix<RealScalar, SizeMinusOne, 1, Options, MaxSizeMinusOne, 1> SubDiagonalType;
typedef Matrix<Scalar, 1, Size, Options, 1, MaxSize> RowVectorType;
typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
typename Diagonal<MatrixType,0>::RealReturnType,
Diagonal<MatrixType,0>
@ -238,7 +240,7 @@ void Tridiagonalization<MatrixType>::matrixQInPlace(MatrixBase<QDerived>* q) con
QDerived& matQ = q->derived();
int n = m_matrix.rows();
matQ = MatrixType::Identity(n,n);
Matrix<Scalar,1,Dynamic> aux(n);
RowVectorType aux(n);
for (int i = n-2; i>=0; i--)
{
matQ.corner(BottomRight,n-i-1,n-i-1)

View File

@ -59,12 +59,19 @@ template<typename _MatrixType> class FullPivLU
{
public:
typedef _MatrixType MatrixType;
enum {
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
Options = MatrixType::Options,
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
};
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> IntRowVectorType;
typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> IntColVectorType;
typedef PermutationMatrix<MatrixType::ColsAtCompileTime> PermutationQType;
typedef PermutationMatrix<MatrixType::RowsAtCompileTime> PermutationPType;
typedef Matrix<int, 1, ColsAtCompileTime, Options, 1, MaxColsAtCompileTime> IntRowVectorType;
typedef Matrix<int, RowsAtCompileTime, 1, Options, MaxRowsAtCompileTime, 1> IntColVectorType;
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationQType;
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationPType;
/**
* \brief Default Constructor.

View File

@ -62,15 +62,18 @@ template<typename _MatrixType> class PartialPivLU
public:
typedef _MatrixType MatrixType;
enum {
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
Options = MatrixType::Options,
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
};
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> PermutationVectorType;
typedef PermutationMatrix<MatrixType::RowsAtCompileTime> PermutationType;
typedef Matrix<int, RowsAtCompileTime, 1, Options, MaxRowsAtCompileTime, 1> PermutationVectorType;
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;
enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN(
MatrixType::MaxColsAtCompileTime,
MatrixType::MaxRowsAtCompileTime)
};
/**
* \brief Default Constructor.

View File

@ -51,16 +51,19 @@ template<typename _MatrixType> class ColPivHouseholderQR
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
Options = MatrixType::Options,
DiagSizeAtCompileTime = EIGEN_SIZE_MIN(ColsAtCompileTime,RowsAtCompileTime)
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
DiagSizeAtCompileTime = EIGEN_SIZE_MIN(ColsAtCompileTime,RowsAtCompileTime),
MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN(MaxColsAtCompileTime,MaxRowsAtCompileTime)
};
typedef typename MatrixType::Scalar Scalar;
typedef typename MatrixType::RealScalar RealScalar;
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime> MatrixQType;
typedef Matrix<Scalar, DiagSizeAtCompileTime, 1> HCoeffsType;
typedef PermutationMatrix<ColsAtCompileTime> PermutationType;
typedef Matrix<int, 1, ColsAtCompileTime> IntRowVectorType;
typedef Matrix<Scalar, 1, ColsAtCompileTime> RowVectorType;
typedef Matrix<RealScalar, 1, ColsAtCompileTime> RealRowVectorType;
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, Options, MaxRowsAtCompileTime, MaxRowsAtCompileTime> MatrixQType;
typedef Matrix<Scalar, DiagSizeAtCompileTime, 1, Options, MaxDiagSizeAtCompileTime, 1> HCoeffsType;
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationType;
typedef Matrix<int, 1, ColsAtCompileTime, Options, 1, MaxColsAtCompileTime> IntRowVectorType;
typedef Matrix<Scalar, 1, ColsAtCompileTime, Options, 1, MaxColsAtCompileTime> RowVectorType;
typedef Matrix<RealScalar, 1, ColsAtCompileTime, Options, 1, MaxColsAtCompileTime> RealRowVectorType;
typedef typename HouseholderSequence<MatrixType,HCoeffsType>::ConjugateReturnType HouseholderSequenceType;
/**

View File

@ -51,17 +51,20 @@ template<typename _MatrixType> class FullPivHouseholderQR
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
Options = MatrixType::Options,
DiagSizeAtCompileTime = EIGEN_SIZE_MIN(ColsAtCompileTime,RowsAtCompileTime)
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
DiagSizeAtCompileTime = EIGEN_SIZE_MIN(ColsAtCompileTime,RowsAtCompileTime),
MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN(MaxColsAtCompileTime,MaxRowsAtCompileTime)
};
typedef typename MatrixType::Scalar Scalar;
typedef typename MatrixType::RealScalar RealScalar;
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime> MatrixQType;
typedef Matrix<Scalar, DiagSizeAtCompileTime, 1> HCoeffsType;
typedef Matrix<int, 1, ColsAtCompileTime> IntRowVectorType;
typedef PermutationMatrix<ColsAtCompileTime> PermutationType;
typedef Matrix<int, RowsAtCompileTime, 1> IntColVectorType;
typedef Matrix<Scalar, 1, ColsAtCompileTime> RowVectorType;
typedef Matrix<Scalar, RowsAtCompileTime, 1> ColVectorType;
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, Options, MaxRowsAtCompileTime, MaxRowsAtCompileTime> MatrixQType;
typedef Matrix<Scalar, DiagSizeAtCompileTime, 1, Options, MaxDiagSizeAtCompileTime, 1> HCoeffsType;
typedef Matrix<int, 1, ColsAtCompileTime, Options, 1, MaxColsAtCompileTime> IntRowVectorType;
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationType;
typedef Matrix<int, RowsAtCompileTime, 1, Options, MaxRowsAtCompileTime, 1> IntColVectorType;
typedef Matrix<Scalar, 1, ColsAtCompileTime, Options, 1, MaxColsAtCompileTime> RowVectorType;
typedef Matrix<Scalar, RowsAtCompileTime, 1, Options, MaxRowsAtCompileTime, 1> ColVectorType;
/** \brief Default Constructor.
*

View File

@ -55,13 +55,16 @@ template<typename _MatrixType> class HouseholderQR
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
Options = MatrixType::Options,
DiagSizeAtCompileTime = EIGEN_SIZE_MIN(ColsAtCompileTime,RowsAtCompileTime)
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
DiagSizeAtCompileTime = EIGEN_SIZE_MIN(ColsAtCompileTime,RowsAtCompileTime),
MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN(MaxColsAtCompileTime,MaxRowsAtCompileTime)
};
typedef typename MatrixType::Scalar Scalar;
typedef typename MatrixType::RealScalar RealScalar;
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, ei_traits<MatrixType>::Flags&RowMajorBit ? RowMajor : ColMajor> MatrixQType;
typedef Matrix<Scalar, DiagSizeAtCompileTime, 1> HCoeffsType;
typedef Matrix<Scalar, 1, ColsAtCompileTime> RowVectorType;
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, ei_traits<MatrixType>::Flags&RowMajorBit ? RowMajor : ColMajor, MaxRowsAtCompileTime, MaxRowsAtCompileTime> MatrixQType;
typedef Matrix<Scalar, DiagSizeAtCompileTime, 1, Options, MaxDiagSizeAtCompileTime, 1> HCoeffsType;
typedef Matrix<Scalar, 1, ColsAtCompileTime, Options, 1, MaxColsAtCompileTime> RowVectorType;
typedef typename HouseholderSequence<MatrixType,HCoeffsType>::ConjugateReturnType HouseholderSequenceType;
/**

View File

@ -52,15 +52,18 @@ template<typename _MatrixType> class SVD
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
PacketSize = ei_packet_traits<Scalar>::size,
AlignmentMask = int(PacketSize)-1,
MinSize = EIGEN_SIZE_MIN(RowsAtCompileTime, ColsAtCompileTime)
MinSize = EIGEN_SIZE_MIN(RowsAtCompileTime, ColsAtCompileTime),
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
MatrixOptions = MatrixType::Options
};
typedef Matrix<Scalar, RowsAtCompileTime, 1> ColVector;
typedef Matrix<Scalar, ColsAtCompileTime, 1> RowVector;
typedef Matrix<Scalar, RowsAtCompileTime, 1, MatrixOptions, MaxRowsAtCompileTime, 1> ColVector;
typedef Matrix<Scalar, ColsAtCompileTime, 1, MatrixOptions, MaxColsAtCompileTime, 1> RowVector;
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime> MatrixUType;
typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime> MatrixVType;
typedef Matrix<Scalar, ColsAtCompileTime, 1> SingularValuesType;
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime> MatrixUType;
typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime> MatrixVType;
typedef Matrix<Scalar, ColsAtCompileTime, 1, MatrixOptions, MaxColsAtCompileTime, 1> SingularValuesType;
/**
* \brief Default Constructor.
@ -195,7 +198,8 @@ SVD<MatrixType>& SVD<MatrixType>::compute(const MatrixType& matrix)
bool convergence = true;
Scalar eps = NumTraits<Scalar>::dummy_precision();
Matrix<Scalar,Dynamic,1> rv1(n);
RowVector rv1(n);
g = scale = anorm = 0;
// Householder reduction to bidiagonal form.
for (i=0; i<n; i++)

View File

@ -33,6 +33,11 @@
#define EIGEN_NO_MALLOC
#include "main.h"
#include <Eigen/Cholesky>
#include <Eigen/Eigenvalues>
#include <Eigen/LU>
#include <Eigen/QR>
#include <Eigen/SVD>
template<typename MatrixType> void nomalloc(const MatrixType& m)
{
@ -73,6 +78,58 @@ template<typename MatrixType> void nomalloc(const MatrixType& m)
}
}
void ctms_decompositions()
{
const int maxSize = 16;
const int size = 12;
typedef Eigen::Matrix<float,
Eigen::Dynamic, Eigen::Dynamic,
Eigen::ColMajor | Eigen::AutoAlign,
maxSize, maxSize> Matrix;
typedef Eigen::Matrix<float,
Eigen::Dynamic, 1,
Eigen::ColMajor | Eigen::AutoAlign,
maxSize, 1> Vector;
typedef Eigen::Matrix<std::complex<float>,
Eigen::Dynamic, Eigen::Dynamic,
Eigen::ColMajor | Eigen::AutoAlign,
maxSize, maxSize> ComplexMatrix;
const Matrix A(Matrix::Random(size, size));
const ComplexMatrix complexA(ComplexMatrix::Random(size, size));
// const Matrix saA = A.adjoint() * A; // NOTE: This product allocates on the stack. The two following lines are a kludgy workaround
Matrix saA(Matrix::Constant(size, size, 1.0));
saA.diagonal().setConstant(2.0);
// Cholesky module
Eigen::LLT<Matrix> LLT; LLT.compute(A);
Eigen::LDLT<Matrix> LDLT; LDLT.compute(A);
// Eigenvalues module
Eigen::HessenbergDecomposition<ComplexMatrix> hessDecomp; hessDecomp.compute(complexA);
Eigen::ComplexSchur<ComplexMatrix> cSchur(size); cSchur.compute(complexA);
Eigen::ComplexEigenSolver<ComplexMatrix> cEigSolver; //cEigSolver.compute(complexA); // NOTE: Commented-out because makes test fail (L135 of ComplexEigenSolver.h has a product that allocates on the stack)
Eigen::EigenSolver<Matrix> eigSolver; eigSolver.compute(A);
Eigen::SelfAdjointEigenSolver<Matrix> saEigSolver(size); saEigSolver.compute(saA);
Eigen::Tridiagonalization<Matrix> tridiag; tridiag.compute(saA);
// LU module
Eigen::PartialPivLU<Matrix> ppLU; ppLU.compute(A);
Eigen::FullPivLU<Matrix> fpLU; fpLU.compute(A);
// QR module
Eigen::HouseholderQR<Matrix> hQR; hQR.compute(A);
Eigen::ColPivHouseholderQR<Matrix> cpQR; cpQR.compute(A);
Eigen::FullPivHouseholderQR<Matrix> fpQR; fpQR.compute(A);
// SVD module
Eigen::JacobiSVD<Matrix> jSVD; jSVD.compute(A);
Eigen::SVD<Matrix> svd; svd.compute(A);
}
void test_nomalloc()
{
// check that our operator new is indeed called:
@ -80,4 +137,8 @@ void test_nomalloc()
CALL_SUBTEST(nomalloc(Matrix<float, 1, 1>()) );
CALL_SUBTEST(nomalloc(Matrix4d()) );
CALL_SUBTEST(nomalloc(Matrix<float,32,32>()) );
// Check decomposition modules with dynamic matrices that have a known compile-time max size (ctms)
CALL_SUBTEST(ctms_decompositions());
}