aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/src/Core/MatrixBase.h2
-rw-r--r--Eigen/src/Core/util/XprHelper.h20
-rw-r--r--Eigen/src/Eigenvalues/ComplexEigenSolver.h36
-rw-r--r--Eigen/src/Geometry/Homogeneous.h14
-rw-r--r--Eigen/src/Householder/BlockHouseholder.h3
-rw-r--r--test/geo_homogeneous.cpp3
-rw-r--r--test/mixingtypes.cpp4
7 files changed, 54 insertions, 28 deletions
diff --git a/Eigen/src/Core/MatrixBase.h b/Eigen/src/Core/MatrixBase.h
index f481c19fb..6674bc65f 100644
--- a/Eigen/src/Core/MatrixBase.h
+++ b/Eigen/src/Core/MatrixBase.h
@@ -129,7 +129,7 @@ template<typename Derived> class MatrixBase
Transpose<Derived>
>::ret AdjointReturnType;
/** \internal Return type of eigenvalues() */
- typedef Matrix<std::complex<RealScalar>, ei_traits<Derived>::ColsAtCompileTime, 1> EigenvaluesReturnType;
+ typedef Matrix<std::complex<RealScalar>, ei_traits<Derived>::ColsAtCompileTime, 1, ColMajor> EigenvaluesReturnType;
/** \internal the return type of identity */
typedef CwiseNullaryOp<ei_scalar_identity_op<Scalar>,Derived> IdentityReturnType;
/** \internal the return type of unit vectors */
diff --git a/Eigen/src/Core/util/XprHelper.h b/Eigen/src/Core/util/XprHelper.h
index 54f15712c..a9cbbfe50 100644
--- a/Eigen/src/Core/util/XprHelper.h
+++ b/Eigen/src/Core/util/XprHelper.h
@@ -91,6 +91,26 @@ template<typename T> struct ei_unpacket_traits
enum {size=1};
};
+template<typename _Scalar, int _Rows, int _Cols,
+ int _Options = AutoAlign |
+ ( (_Rows==1 && _Cols!=1) ? RowMajor
+ : (_Cols==1 && _Rows!=1) ? ColMajor
+ : EIGEN_DEFAULT_MATRIX_STORAGE_ORDER_OPTION ),
+ int _MaxRows = _Rows,
+ int _MaxCols = _Cols
+> class ei_make_proper_matrix_type
+{
+ enum {
+ IsColVector = _Cols==1 && _Rows!=1,
+ IsRowVector = _Rows==1 && _Cols!=1,
+ Options = IsColVector ? (_Options | ColMajor) & ~RowMajor
+ : IsRowVector ? (_Options | RowMajor) & ~ColMajor
+ : _Options
+ };
+ public:
+ typedef Matrix<_Scalar, _Rows, _Cols, Options, _MaxRows, _MaxCols> type;
+};
+
template<typename Scalar, int Rows, int Cols, int Options, int MaxRows, int MaxCols>
class ei_compute_matrix_flags
{
diff --git a/Eigen/src/Eigenvalues/ComplexEigenSolver.h b/Eigen/src/Eigenvalues/ComplexEigenSolver.h
index a164aaae6..f9e53c77a 100644
--- a/Eigen/src/Eigenvalues/ComplexEigenSolver.h
+++ b/Eigen/src/Eigenvalues/ComplexEigenSolver.h
@@ -76,7 +76,7 @@ template<typename _MatrixType> class ComplexEigenSolver
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef typename MatrixType::Index Index;
- /** \brief Complex scalar type for #MatrixType.
+ /** \brief Complex scalar type for #MatrixType.
*
* This is \c std::complex<Scalar> if #Scalar is real (e.g.,
* \c float or \c double) and just \c Scalar if #Scalar is
@@ -84,16 +84,16 @@ template<typename _MatrixType> class ComplexEigenSolver
*/
typedef std::complex<RealScalar> ComplexScalar;
- /** \brief Type for vector of eigenvalues as returned by eigenvalues().
+ /** \brief Type for vector of eigenvalues as returned by eigenvalues().
*
* This is a column vector with entries of type #ComplexScalar.
* The length of the vector is the size of #MatrixType.
*/
- typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options, MaxColsAtCompileTime, 1> EigenvalueType;
+ typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options&(~RowMajor), MaxColsAtCompileTime, 1> EigenvalueType;
- /** \brief Type for matrix of eigenvectors as returned by eigenvectors().
+ /** \brief Type for matrix of eigenvectors as returned by eigenvectors().
*
- * This is a square matrix with entries of type #ComplexScalar.
+ * This is a square matrix with entries of type #ComplexScalar.
* The size is the same as the size of #MatrixType.
*/
typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, ColsAtCompileTime> EigenvectorType;
@@ -111,7 +111,7 @@ template<typename _MatrixType> class ComplexEigenSolver
m_eigenvectorsOk(false),
m_matX()
{}
-
+
/** \brief Default Constructor with memory preallocation
*
* Like the default constructor but with preallocation of the internal data
@@ -127,12 +127,12 @@ template<typename _MatrixType> class ComplexEigenSolver
m_matX(size, size)
{}
- /** \brief Constructor; computes eigendecomposition of given matrix.
- *
+ /** \brief Constructor; computes eigendecomposition of given matrix.
+ *
* \param[in] matrix Square matrix whose eigendecomposition is to be computed.
* \param[in] computeEigenvectors If true, both the eigenvectors and the
* eigenvalues are computed; if false, only the eigenvalues are
- * computed.
+ * computed.
*
* This constructor calls compute() to compute the eigendecomposition.
*/
@@ -147,14 +147,14 @@ template<typename _MatrixType> class ComplexEigenSolver
compute(matrix, computeEigenvectors);
}
- /** \brief Returns the eigenvectors of given matrix.
+ /** \brief Returns the eigenvectors of given matrix.
*
* \returns A const reference to the matrix whose columns are the eigenvectors.
*
* \pre Either the constructor
* ComplexEigenSolver(const MatrixType& matrix, bool) or the member
* function compute(const MatrixType& matrix, bool) has been called before
- * to compute the eigendecomposition of a matrix, and
+ * to compute the eigendecomposition of a matrix, and
* \p computeEigenvectors was set to true (the default).
*
* This function returns a matrix whose columns are the eigenvectors. Column
@@ -174,7 +174,7 @@ template<typename _MatrixType> class ComplexEigenSolver
return m_eivec;
}
- /** \brief Returns the eigenvalues of given matrix.
+ /** \brief Returns the eigenvalues of given matrix.
*
* \returns A const reference to the column vector containing the eigenvalues.
*
@@ -197,16 +197,16 @@ template<typename _MatrixType> class ComplexEigenSolver
return m_eivalues;
}
- /** \brief Computes eigendecomposition of given matrix.
- *
+ /** \brief Computes eigendecomposition of given matrix.
+ *
* \param[in] matrix Square matrix whose eigendecomposition is to be computed.
* \param[in] computeEigenvectors If true, both the eigenvectors and the
* eigenvalues are computed; if false, only the eigenvalues are
- * computed.
+ * computed.
* \returns Reference to \c *this
*
* This function computes the eigenvalues of the complex matrix \p matrix.
- * The eigenvalues() function can be used to retrieve them. If
+ * The eigenvalues() function can be used to retrieve them. If
* \p computeEigenvectors is true, then the eigenvectors are also computed
* and can be retrieved by calling eigenvectors().
*
@@ -257,7 +257,7 @@ ComplexEigenSolver<MatrixType>& ComplexEigenSolver<MatrixType>::compute(const Ma
// The eigenvalues are on the diagonal of T.
m_schur.compute(matrix, computeEigenvectors);
- if(m_schur.info() == Success)
+ if(m_schur.info() == Success)
{
m_eivalues = m_schur.matrixT().diagonal();
if(computeEigenvectors)
@@ -291,7 +291,7 @@ void ComplexEigenSolver<MatrixType>::doComputeEigenvectors(RealScalar matrixnorm
ComplexScalar z = m_schur.matrixT().coeff(i,i) - m_schur.matrixT().coeff(k,k);
if(z==ComplexScalar(0))
{
- // If the i-th and k-th eigenvalue are equal, then z equals 0.
+ // If the i-th and k-th eigenvalue are equal, then z equals 0.
// Use a small value instead, to prevent division by zero.
ei_real_ref(z) = NumTraits<RealScalar>::epsilon() * matrixnorm;
}
diff --git a/Eigen/src/Geometry/Homogeneous.h b/Eigen/src/Geometry/Homogeneous.h
index 3077f0921..fc6134fe3 100644
--- a/Eigen/src/Geometry/Homogeneous.h
+++ b/Eigen/src/Geometry/Homogeneous.h
@@ -55,7 +55,10 @@ struct ei_traits<Homogeneous<MatrixType,Direction> >
ColsAtCompileTime = Direction==Horizontal ? ColsPlusOne : MatrixType::ColsAtCompileTime,
MaxRowsAtCompileTime = RowsAtCompileTime,
MaxColsAtCompileTime = ColsAtCompileTime,
- Flags = _MatrixTypeNested::Flags & HereditaryBits,
+ TmpFlags = _MatrixTypeNested::Flags & HereditaryBits,
+ Flags = ColsAtCompileTime==1 ? (TmpFlags & ~RowMajorBit)
+ : RowsAtCompileTime==1 ? (TmpFlags | RowMajorBit)
+ : TmpFlags,
CoeffReadCost = _MatrixTypeNested::CoeffReadCost
};
};
@@ -210,12 +213,13 @@ VectorwiseOp<ExpressionType,Direction>::hnormalized() const
template<typename MatrixType,typename Lhs>
struct ei_traits<ei_homogeneous_left_product_impl<Homogeneous<MatrixType,Vertical>,Lhs> >
{
- typedef Matrix<typename ei_traits<MatrixType>::Scalar,
+ typedef typename ei_make_proper_matrix_type<
+ typename ei_traits<MatrixType>::Scalar,
Lhs::RowsAtCompileTime,
MatrixType::ColsAtCompileTime,
MatrixType::PlainObject::Options,
Lhs::MaxRowsAtCompileTime,
- MatrixType::MaxColsAtCompileTime> ReturnType;
+ MatrixType::MaxColsAtCompileTime>::type ReturnType;
};
template<typename MatrixType,typename Lhs>
@@ -249,12 +253,12 @@ struct ei_homogeneous_left_product_impl<Homogeneous<MatrixType,Vertical>,Lhs>
template<typename MatrixType,typename Rhs>
struct ei_traits<ei_homogeneous_right_product_impl<Homogeneous<MatrixType,Horizontal>,Rhs> >
{
- typedef Matrix<typename ei_traits<MatrixType>::Scalar,
+ typedef typename ei_make_proper_matrix_type<typename ei_traits<MatrixType>::Scalar,
MatrixType::RowsAtCompileTime,
Rhs::ColsAtCompileTime,
MatrixType::PlainObject::Options,
MatrixType::MaxRowsAtCompileTime,
- Rhs::MaxColsAtCompileTime> ReturnType;
+ Rhs::MaxColsAtCompileTime>::type ReturnType;
};
template<typename MatrixType,typename Rhs>
diff --git a/Eigen/src/Householder/BlockHouseholder.h b/Eigen/src/Householder/BlockHouseholder.h
index c17ca10e0..60da512e1 100644
--- a/Eigen/src/Householder/BlockHouseholder.h
+++ b/Eigen/src/Householder/BlockHouseholder.h
@@ -66,7 +66,8 @@ void ei_apply_block_householder_on_the_left(MatrixType& mat, const VectorsType&
const TriangularView<VectorsType, UnitLower>& V(vectors);
// A -= V T V^* A
- Matrix<typename MatrixType::Scalar,Dynamic,Dynamic> tmp = V.adjoint() * mat;
+ Matrix<typename MatrixType::Scalar,VectorsType::ColsAtCompileTime,MatrixType::ColsAtCompileTime,0,
+ VectorsType::MaxColsAtCompileTime,MatrixType::MaxColsAtCompileTime> tmp = V.adjoint() * mat;
// FIXME add .noalias() once the triangular product can work inplace
tmp = T.template triangularView<Upper>().adjoint() * tmp;
mat.noalias() -= V * tmp;
diff --git a/test/geo_homogeneous.cpp b/test/geo_homogeneous.cpp
index 781913553..3162af831 100644
--- a/test/geo_homogeneous.cpp
+++ b/test/geo_homogeneous.cpp
@@ -32,7 +32,7 @@ template<typename Scalar,int Size> void homogeneous(void)
*/
typedef Matrix<Scalar,Size,Size> MatrixType;
- typedef Matrix<Scalar,Size,1> VectorType;
+ typedef Matrix<Scalar,Size,1, ColMajor> VectorType;
typedef Matrix<Scalar,Size+1,Size> HMatrixType;
typedef Matrix<Scalar,Size+1,1> HVectorType;
@@ -80,6 +80,7 @@ template<typename Scalar,int Size> void homogeneous(void)
VERIFY_IS_APPROX((v0.transpose().rowwise().homogeneous().eval()) * t2,
v0.transpose().rowwise().homogeneous() * t2);
+ m0.transpose().rowwise().homogeneous().eval();
VERIFY_IS_APPROX((m0.transpose().rowwise().homogeneous().eval()) * t2,
m0.transpose().rowwise().homogeneous() * t2);
diff --git a/test/mixingtypes.cpp b/test/mixingtypes.cpp
index 0465f3059..3b537e89d 100644
--- a/test/mixingtypes.cpp
+++ b/test/mixingtypes.cpp
@@ -136,12 +136,12 @@ void mixingtypes_large(int size)
VERIFY_RAISES_ASSERT(mcf*vf);
// VERIFY_RAISES_ASSERT(mcf *= mf); // does not even compile
// VERIFY_RAISES_ASSERT(vcd = md*vcd); // does not even compile (cannot convert complex to double)
- VERIFY_RAISES_ASSERT(vcf = mcf*vf);
+// VERIFY_RAISES_ASSERT(vcf = mcf*vf);
// VERIFY_RAISES_ASSERT(mf*md); // does not even compile
// VERIFY_RAISES_ASSERT(mcf*mcd); // does not even compile
// VERIFY_RAISES_ASSERT(mcf*vcd); // does not even compile
- VERIFY_RAISES_ASSERT(vcf = mf*vf);
+// VERIFY_RAISES_ASSERT(vcf = mf*vf);
}
template<int SizeAtCompileType> void mixingtypes_small()