aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Core
diff options
context:
space:
mode:
authorGravatar Hauke Heibel <hauke.heibel@gmail.com>2010-02-20 15:53:57 +0100
committerGravatar Hauke Heibel <hauke.heibel@gmail.com>2010-02-20 15:53:57 +0100
commitabc8c010807f0706b80bc0ef13303b6485473a57 (patch)
tree1767172943f08673a1df661273b97b6bbd48c546 /Eigen/src/Core
parentf0c8dcf1e2f01fb200a8e48463d9f73c77bc1436 (diff)
Renamed PlainMatrixType to PlainObject (Array != Matrix).
Renamed ReturnByValue::ReturnMatrixType ReturnByValue::ReturnType (again, Array != Matrix).
Diffstat (limited to 'Eigen/src/Core')
-rw-r--r--Eigen/src/Core/DenseStorageBase.h10
-rw-r--r--Eigen/src/Core/Dot.h2
-rw-r--r--Eigen/src/Core/EigenBase.h6
-rw-r--r--Eigen/src/Core/Flagged.h2
-rw-r--r--Eigen/src/Core/Matrix.h2
-rw-r--r--Eigen/src/Core/MatrixBase.h31
-rw-r--r--Eigen/src/Core/PermutationMatrix.h2
-rw-r--r--Eigen/src/Core/ProductBase.h14
-rw-r--r--Eigen/src/Core/ReturnByValue.h16
-rw-r--r--Eigen/src/Core/SelfAdjointView.h6
-rw-r--r--Eigen/src/Core/SelfCwiseBinaryOp.h8
-rw-r--r--Eigen/src/Core/TriangularMatrix.h2
-rw-r--r--Eigen/src/Core/products/CoeffBasedProduct.h14
-rw-r--r--Eigen/src/Core/util/BlasUtil.h4
-rw-r--r--Eigen/src/Core/util/XprHelper.h6
15 files changed, 62 insertions, 63 deletions
diff --git a/Eigen/src/Core/DenseStorageBase.h b/Eigen/src/Core/DenseStorageBase.h
index 530ddcd07..6245b8007 100644
--- a/Eigen/src/Core/DenseStorageBase.h
+++ b/Eigen/src/Core/DenseStorageBase.h
@@ -44,7 +44,7 @@ class DenseStorageBase : public _Base<Derived>
public:
enum { Options = _Options };
typedef _Base<Derived> Base;
- typedef typename Base::PlainMatrixType PlainMatrixType;
+ typedef typename Base::PlainObject PlainObject;
typedef typename Base::Scalar Scalar;
typedef typename Base::PacketScalar PacketScalar;
using Base::RowsAtCompileTime;
@@ -544,7 +544,7 @@ struct ei_conservative_resize_like_impl
{
if (_this.rows() == rows && _this.cols() == cols) return;
EIGEN_STATIC_ASSERT_DYNAMIC_SIZE(Derived)
- typename Derived::PlainMatrixType tmp(rows,cols);
+ typename Derived::PlainObject tmp(rows,cols);
const int common_rows = std::min(rows, _this.rows());
const int common_cols = std::min(cols, _this.cols());
tmp.block(0,0,common_rows,common_cols) = _this.block(0,0,common_rows,common_cols);
@@ -563,7 +563,7 @@ struct ei_conservative_resize_like_impl
EIGEN_STATIC_ASSERT_DYNAMIC_SIZE(Derived)
EIGEN_STATIC_ASSERT_DYNAMIC_SIZE(OtherDerived)
- typename Derived::PlainMatrixType tmp(other);
+ typename Derived::PlainObject tmp(other);
const int common_rows = std::min(tmp.rows(), _this.rows());
const int common_cols = std::min(tmp.cols(), _this.cols());
tmp.block(0,0,common_rows,common_cols) = _this.block(0,0,common_rows,common_cols);
@@ -577,7 +577,7 @@ struct ei_conservative_resize_like_impl<Derived,OtherDerived,true>
static void run(DenseBase<Derived>& _this, int size)
{
if (_this.size() == size) return;
- typename Derived::PlainMatrixType tmp(size);
+ typename Derived::PlainObject tmp(size);
const int common_size = std::min<int>(_this.size(),size);
tmp.segment(0,common_size) = _this.segment(0,common_size);
_this.derived().swap(tmp);
@@ -588,7 +588,7 @@ struct ei_conservative_resize_like_impl<Derived,OtherDerived,true>
if (_this.rows() == other.rows() && _this.cols() == other.cols()) return;
// segment(...) will check whether Derived/OtherDerived are vectors!
- typename Derived::PlainMatrixType tmp(other);
+ typename Derived::PlainObject tmp(other);
const int common_size = std::min<int>(_this.size(),tmp.size());
tmp.segment(0,common_size) = _this.segment(0,common_size);
_this.derived().swap(tmp);
diff --git a/Eigen/src/Core/Dot.h b/Eigen/src/Core/Dot.h
index f0c520b1f..201bd23ca 100644
--- a/Eigen/src/Core/Dot.h
+++ b/Eigen/src/Core/Dot.h
@@ -299,7 +299,7 @@ inline typename NumTraits<typename ei_traits<Derived>::Scalar>::Real MatrixBase<
* \sa norm(), normalize()
*/
template<typename Derived>
-inline const typename MatrixBase<Derived>::PlainMatrixType
+inline const typename MatrixBase<Derived>::PlainObject
MatrixBase<Derived>::normalized() const
{
typedef typename ei_nested<Derived>::type Nested;
diff --git a/Eigen/src/Core/EigenBase.h b/Eigen/src/Core/EigenBase.h
index 8b302b663..cf1ce4376 100644
--- a/Eigen/src/Core/EigenBase.h
+++ b/Eigen/src/Core/EigenBase.h
@@ -37,7 +37,7 @@
*/
template<typename Derived> struct EigenBase
{
-// typedef typename ei_plain_matrix_type<Derived>::type PlainMatrixType;
+// typedef typename ei_plain_matrix_type<Derived>::type PlainObject;
/** \returns a reference to the derived object */
Derived& derived() { return *static_cast<Derived*>(this); }
@@ -61,7 +61,7 @@ template<typename Derived> struct EigenBase
{
// This is the default implementation,
// derived class can reimplement it in a more optimized way.
- typename Dest::PlainMatrixType res(rows(),cols());
+ typename Dest::PlainObject res(rows(),cols());
evalTo(res);
dst += res;
}
@@ -71,7 +71,7 @@ template<typename Derived> struct EigenBase
{
// This is the default implementation,
// derived class can reimplement it in a more optimized way.
- typename Dest::PlainMatrixType res(rows(),cols());
+ typename Dest::PlainObject res(rows(),cols());
evalTo(res);
dst -= res;
}
diff --git a/Eigen/src/Core/Flagged.h b/Eigen/src/Core/Flagged.h
index 7f42a1e73..0044fe7cb 100644
--- a/Eigen/src/Core/Flagged.h
+++ b/Eigen/src/Core/Flagged.h
@@ -109,7 +109,7 @@ template<typename ExpressionType, unsigned int Added, unsigned int Removed> clas
const ExpressionType& _expression() const { return m_matrix; }
template<typename OtherDerived>
- typename ExpressionType::PlainMatrixType solveTriangular(const MatrixBase<OtherDerived>& other) const;
+ typename ExpressionType::PlainObject solveTriangular(const MatrixBase<OtherDerived>& other) const;
template<typename OtherDerived>
void solveTriangularInPlace(const MatrixBase<OtherDerived>& other) const;
diff --git a/Eigen/src/Core/Matrix.h b/Eigen/src/Core/Matrix.h
index 1c43340a6..44a0ef7d1 100644
--- a/Eigen/src/Core/Matrix.h
+++ b/Eigen/src/Core/Matrix.h
@@ -139,7 +139,7 @@ class Matrix
EIGEN_DENSE_PUBLIC_INTERFACE(Matrix)
- typedef typename Base::PlainMatrixType PlainMatrixType;
+ typedef typename Base::PlainObject PlainObject;
enum { NeedsToAlign = (!(Options&DontAlign))
&& SizeAtCompileTime!=Dynamic && ((sizeof(Scalar)*SizeAtCompileTime)%16)==0 };
diff --git a/Eigen/src/Core/MatrixBase.h b/Eigen/src/Core/MatrixBase.h
index 122a2271b..9c62163ba 100644
--- a/Eigen/src/Core/MatrixBase.h
+++ b/Eigen/src/Core/MatrixBase.h
@@ -121,7 +121,7 @@ template<typename Derived> class MatrixBase
*
* This is not necessarily exactly the return type of eval(). In the case of plain matrices,
* the return type of eval() is a const reference to a matrix, not a matrix! It is however guaranteed
- * that the return type of eval() is either PlainMatrixType or const PlainMatrixType&.
+ * that the return type of eval() is either PlainObject or const PlainObject&.
*/
typedef Matrix<typename ei_traits<Derived>::Scalar,
ei_traits<Derived>::RowsAtCompileTime,
@@ -129,8 +129,7 @@ template<typename Derived> class MatrixBase
AutoAlign | (ei_traits<Derived>::Flags&RowMajorBit ? RowMajor : ColMajor),
ei_traits<Derived>::MaxRowsAtCompileTime,
ei_traits<Derived>::MaxColsAtCompileTime
- > PlainMatrixType;
- // typedef typename ei_plain_matrix_type<Derived>::type PlainMatrixType;
+ > PlainObject;
#ifndef EIGEN_PARSED_BY_DOXYGEN
/** \internal Represents a matrix with all coefficients equal to one another*/
@@ -212,7 +211,7 @@ template<typename Derived> class MatrixBase
RealScalar stableNorm() const;
RealScalar blueNorm() const;
RealScalar hypotNorm() const;
- const PlainMatrixType normalized() const;
+ const PlainObject normalized() const;
void normalize();
const AdjointReturnType adjoint() const;
@@ -301,9 +300,9 @@ template<typename Derived> class MatrixBase
/////////// LU module ///////////
- const FullPivLU<PlainMatrixType> fullPivLu() const;
- const PartialPivLU<PlainMatrixType> partialPivLu() const;
- const PartialPivLU<PlainMatrixType> lu() const;
+ const FullPivLU<PlainObject> fullPivLu() const;
+ const PartialPivLU<PlainObject> partialPivLu() const;
+ const PartialPivLU<PlainObject> lu() const;
const ei_inverse_impl<Derived> inverse() const;
template<typename ResultType>
void computeInverseAndDetWithCheck(
@@ -322,29 +321,29 @@ template<typename Derived> class MatrixBase
/////////// Cholesky module ///////////
- const LLT<PlainMatrixType> llt() const;
- const LDLT<PlainMatrixType> ldlt() const;
+ const LLT<PlainObject> llt() const;
+ const LDLT<PlainObject> ldlt() const;
/////////// QR module ///////////
- const HouseholderQR<PlainMatrixType> householderQr() const;
- const ColPivHouseholderQR<PlainMatrixType> colPivHouseholderQr() const;
- const FullPivHouseholderQR<PlainMatrixType> fullPivHouseholderQr() const;
+ const HouseholderQR<PlainObject> householderQr() const;
+ const ColPivHouseholderQR<PlainObject> colPivHouseholderQr() const;
+ const FullPivHouseholderQR<PlainObject> fullPivHouseholderQr() const;
EigenvaluesReturnType eigenvalues() const;
RealScalar operatorNorm() const;
/////////// SVD module ///////////
- SVD<PlainMatrixType> svd() const;
+ SVD<PlainObject> svd() const;
/////////// Geometry module ///////////
template<typename OtherDerived>
- PlainMatrixType cross(const MatrixBase<OtherDerived>& other) const;
+ PlainObject cross(const MatrixBase<OtherDerived>& other) const;
template<typename OtherDerived>
- PlainMatrixType cross3(const MatrixBase<OtherDerived>& other) const;
- PlainMatrixType unitOrthogonal(void) const;
+ PlainObject cross3(const MatrixBase<OtherDerived>& other) const;
+ PlainObject unitOrthogonal(void) const;
Matrix<Scalar,3,1> eulerAngles(int a0, int a1, int a2) const;
const ScalarMultipleReturnType operator*(const UniformScaling<Scalar>& s) const;
enum {
diff --git a/Eigen/src/Core/PermutationMatrix.h b/Eigen/src/Core/PermutationMatrix.h
index 2d97c9c38..fcd2e46cc 100644
--- a/Eigen/src/Core/PermutationMatrix.h
+++ b/Eigen/src/Core/PermutationMatrix.h
@@ -280,7 +280,7 @@ operator*(const PermutationMatrix<SizeAtCompileTime, MaxSizeAtCompileTime> &perm
template<typename PermutationType, typename MatrixType, int Side>
struct ei_traits<ei_permut_matrix_product_retval<PermutationType, MatrixType, Side> >
{
- typedef typename MatrixType::PlainMatrixType ReturnMatrixType;
+ typedef typename MatrixType::PlainObject ReturnType;
};
template<typename PermutationType, typename MatrixType, int Side>
diff --git a/Eigen/src/Core/ProductBase.h b/Eigen/src/Core/ProductBase.h
index 481e7c760..789aecfb6 100644
--- a/Eigen/src/Core/ProductBase.h
+++ b/Eigen/src/Core/ProductBase.h
@@ -88,7 +88,7 @@ class ProductBase : public MatrixBase<Derived>
public:
- typedef typename Base::PlainMatrixType PlainMatrixType;
+ typedef typename Base::PlainObject PlainObject;
ProductBase(const Lhs& lhs, const Rhs& rhs)
: m_lhs(lhs), m_rhs(rhs)
@@ -116,8 +116,8 @@ class ProductBase : public MatrixBase<Derived>
const _LhsNested& lhs() const { return m_lhs; }
const _RhsNested& rhs() const { return m_rhs; }
- // Implicit convertion to the nested type (trigger the evaluation of the product)
- operator const PlainMatrixType& () const
+ // Implicit conversion to the nested type (trigger the evaluation of the product)
+ operator const PlainObject& () const
{
m_result.resize(m_lhs.rows(), m_rhs.cols());
this->evalTo(m_result);
@@ -139,7 +139,7 @@ class ProductBase : public MatrixBase<Derived>
const LhsNested m_lhs;
const RhsNested m_rhs;
- mutable PlainMatrixType m_result;
+ mutable PlainObject m_result;
private:
@@ -152,10 +152,10 @@ class ProductBase : public MatrixBase<Derived>
// here we need to overload the nested rule for products
// such that the nested type is a const reference to a plain matrix
-template<typename Lhs, typename Rhs, int Mode, int N, typename PlainMatrixType>
-struct ei_nested<GeneralProduct<Lhs,Rhs,Mode>, N, PlainMatrixType>
+template<typename Lhs, typename Rhs, int Mode, int N, typename PlainObject>
+struct ei_nested<GeneralProduct<Lhs,Rhs,Mode>, N, PlainObject>
{
- typedef PlainMatrixType const& type;
+ typedef PlainObject const& type;
};
template<typename NestedProduct>
diff --git a/Eigen/src/Core/ReturnByValue.h b/Eigen/src/Core/ReturnByValue.h
index 920269365..d375f0b5c 100644
--- a/Eigen/src/Core/ReturnByValue.h
+++ b/Eigen/src/Core/ReturnByValue.h
@@ -31,7 +31,7 @@
*/
template<typename Derived>
struct ei_traits<ReturnByValue<Derived> >
- : public ei_traits<typename ei_traits<Derived>::ReturnMatrixType>
+ : public ei_traits<typename ei_traits<Derived>::ReturnType>
{
enum {
// FIXME had to remove the DirectAccessBit for usage like
@@ -42,7 +42,7 @@ struct ei_traits<ReturnByValue<Derived> >
// The fact that I had to do that shows that when doing xpr.block() with a non-direct-access xpr,
// even if xpr has the EvalBeforeNestingBit, the block() doesn't use direct access on the evaluated
// xpr.
- Flags = (ei_traits<typename ei_traits<Derived>::ReturnMatrixType>::Flags
+ Flags = (ei_traits<typename ei_traits<Derived>::ReturnType>::Flags
| EvalBeforeNestingBit) & ~DirectAccessBit
};
};
@@ -51,18 +51,18 @@ struct ei_traits<ReturnByValue<Derived> >
* So the only way that nesting it in an expression can work, is by evaluating it into a plain matrix.
* So ei_nested always gives the plain return matrix type.
*/
-template<typename Derived,int n,typename PlainMatrixType>
-struct ei_nested<ReturnByValue<Derived>, n, PlainMatrixType>
+template<typename Derived,int n,typename PlainObject>
+struct ei_nested<ReturnByValue<Derived>, n, PlainObject>
{
- typedef typename ei_traits<Derived>::ReturnMatrixType type;
+ typedef typename ei_traits<Derived>::ReturnType type;
};
template<typename Derived> class ReturnByValue
- : public ei_traits<Derived>::ReturnMatrixType::template MakeBase<ReturnByValue<Derived> >::Type
+ : public ei_traits<Derived>::ReturnType::template MakeBase<ReturnByValue<Derived> >::Type
{
public:
- typedef typename ei_traits<Derived>::ReturnMatrixType ReturnMatrixType;
- typedef typename ReturnMatrixType::template MakeBase<ReturnByValue<Derived> >::Type Base;
+ typedef typename ei_traits<Derived>::ReturnType ReturnType;
+ typedef typename ReturnType::template MakeBase<ReturnByValue<Derived> >::Type Base;
EIGEN_DENSE_PUBLIC_INTERFACE(ReturnByValue)
template<typename Dest>
diff --git a/Eigen/src/Core/SelfAdjointView.h b/Eigen/src/Core/SelfAdjointView.h
index 6d01ee495..c3c5b17ff 100644
--- a/Eigen/src/Core/SelfAdjointView.h
+++ b/Eigen/src/Core/SelfAdjointView.h
@@ -68,7 +68,7 @@ template<typename MatrixType, unsigned int UpLo> class SelfAdjointView
enum {
Mode = ei_traits<SelfAdjointView>::Mode
};
- typedef typename MatrixType::PlainMatrixType PlainMatrixType;
+ typedef typename MatrixType::PlainObject PlainObject;
inline SelfAdjointView(const MatrixType& matrix) : m_matrix(matrix)
{ ei_assert(ei_are_flags_consistent<Mode>::ret); }
@@ -146,8 +146,8 @@ template<typename MatrixType, unsigned int UpLo> class SelfAdjointView
/////////// Cholesky module ///////////
- const LLT<PlainMatrixType, UpLo> llt() const;
- const LDLT<PlainMatrixType> ldlt() const;
+ const LLT<PlainObject, UpLo> llt() const;
+ const LDLT<PlainObject> ldlt() const;
protected:
const typename MatrixType::Nested m_matrix;
diff --git a/Eigen/src/Core/SelfCwiseBinaryOp.h b/Eigen/src/Core/SelfCwiseBinaryOp.h
index 7ae2e82a4..d2690b66b 100644
--- a/Eigen/src/Core/SelfCwiseBinaryOp.h
+++ b/Eigen/src/Core/SelfCwiseBinaryOp.h
@@ -124,8 +124,8 @@ template<typename Derived>
inline Derived& DenseBase<Derived>::operator*=(const Scalar& other)
{
SelfCwiseBinaryOp<ei_scalar_product_op<Scalar>, Derived> tmp(derived());
- typedef typename Derived::PlainMatrixType PlainMatrixType;
- tmp = PlainMatrixType::Constant(rows(),cols(),other);
+ typedef typename Derived::PlainObject PlainObject;
+ tmp = PlainObject::Constant(rows(),cols(),other);
return derived();
}
@@ -133,8 +133,8 @@ template<typename Derived>
inline Derived& DenseBase<Derived>::operator/=(const Scalar& other)
{
SelfCwiseBinaryOp<typename ei_meta_if<NumTraits<Scalar>::HasFloatingPoint,ei_scalar_product_op<Scalar>,ei_scalar_quotient_op<Scalar> >::ret, Derived> tmp(derived());
- typedef typename Derived::PlainMatrixType PlainMatrixType;
- tmp = PlainMatrixType::Constant(rows(),cols(), NumTraits<Scalar>::HasFloatingPoint ? Scalar(1)/other : other);
+ typedef typename Derived::PlainObject PlainObject;
+ tmp = PlainObject::Constant(rows(),cols(), NumTraits<Scalar>::HasFloatingPoint ? Scalar(1)/other : other);
return derived();
}
diff --git a/Eigen/src/Core/TriangularMatrix.h b/Eigen/src/Core/TriangularMatrix.h
index 859f10298..7d978b800 100644
--- a/Eigen/src/Core/TriangularMatrix.h
+++ b/Eigen/src/Core/TriangularMatrix.h
@@ -148,7 +148,7 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularView
typedef TriangularBase<TriangularView> Base;
typedef typename ei_traits<TriangularView>::Scalar Scalar;
typedef _MatrixType MatrixType;
- typedef typename MatrixType::PlainMatrixType DenseMatrixType;
+ typedef typename MatrixType::PlainObject DenseMatrixType;
typedef typename MatrixType::Nested MatrixTypeNested;
typedef typename ei_cleantype<MatrixTypeNested>::type _MatrixTypeNested;
diff --git a/Eigen/src/Core/products/CoeffBasedProduct.h b/Eigen/src/Core/products/CoeffBasedProduct.h
index f030d59b5..3343b1875 100644
--- a/Eigen/src/Core/products/CoeffBasedProduct.h
+++ b/Eigen/src/Core/products/CoeffBasedProduct.h
@@ -109,7 +109,7 @@ class CoeffBasedProduct
typedef MatrixBase<CoeffBasedProduct> Base;
EIGEN_DENSE_PUBLIC_INTERFACE(CoeffBasedProduct)
- typedef typename Base::PlainMatrixType PlainMatrixType;
+ typedef typename Base::PlainObject PlainObject;
private:
@@ -181,8 +181,8 @@ class CoeffBasedProduct
return res;
}
- // Implicit convertion to the nested type (trigger the evaluation of the product)
- operator const PlainMatrixType& () const
+ // Implicit conversion to the nested type (trigger the evaluation of the product)
+ operator const PlainObject& () const
{
m_result.lazyAssign(*this);
return m_result;
@@ -205,15 +205,15 @@ class CoeffBasedProduct
const LhsNested m_lhs;
const RhsNested m_rhs;
- mutable PlainMatrixType m_result;
+ mutable PlainObject m_result;
};
// here we need to overload the nested rule for products
// such that the nested type is a const reference to a plain matrix
-template<typename Lhs, typename Rhs, int N, typename PlainMatrixType>
-struct ei_nested<CoeffBasedProduct<Lhs,Rhs,EvalBeforeNestingBit|EvalBeforeAssigningBit>, N, PlainMatrixType>
+template<typename Lhs, typename Rhs, int N, typename PlainObject>
+struct ei_nested<CoeffBasedProduct<Lhs,Rhs,EvalBeforeNestingBit|EvalBeforeAssigningBit>, N, PlainObject>
{
- typedef PlainMatrixType const& type;
+ typedef PlainObject const& type;
};
/***************************************************************************
diff --git a/Eigen/src/Core/util/BlasUtil.h b/Eigen/src/Core/util/BlasUtil.h
index 3777464dc..2ca463d5d 100644
--- a/Eigen/src/Core/util/BlasUtil.h
+++ b/Eigen/src/Core/util/BlasUtil.h
@@ -166,7 +166,7 @@ template<typename XprType> struct ei_blas_traits
};
typedef typename ei_meta_if<int(ActualAccess)==HasDirectAccess,
ExtractType,
- typename _ExtractType::PlainMatrixType
+ typename _ExtractType::PlainObject
>::ret DirectLinearAccessType;
static inline ExtractType extract(const XprType& x) { return x; }
static inline Scalar extractScalarFactor(const XprType&) { return Scalar(1); }
@@ -227,7 +227,7 @@ struct ei_blas_traits<Transpose<NestedXpr> >
typedef Transpose<typename Base::_ExtractType> _ExtractType;
typedef typename ei_meta_if<int(Base::ActualAccess)==HasDirectAccess,
ExtractType,
- typename ExtractType::PlainMatrixType
+ typename ExtractType::PlainObject
>::ret DirectLinearAccessType;
enum {
IsTransposed = Base::IsTransposed ? 0 : 1
diff --git a/Eigen/src/Core/util/XprHelper.h b/Eigen/src/Core/util/XprHelper.h
index 4884557e5..a54ddd155 100644
--- a/Eigen/src/Core/util/XprHelper.h
+++ b/Eigen/src/Core/util/XprHelper.h
@@ -147,7 +147,7 @@ template<typename T, typename StorageType = typename ei_traits<T>::StorageType>
template<typename T> struct ei_eval<T,Dense>
{
typedef typename ei_plain_matrix_type<T>::type type;
-// typedef typename T::PlainMatrixType type;
+// typedef typename T::PlainObject type;
// typedef T::Matrix<typename ei_traits<T>::Scalar,
// ei_traits<T>::RowsAtCompileTime,
// ei_traits<T>::ColsAtCompileTime,
@@ -256,7 +256,7 @@ struct ei_ref_selector
* const Matrix3d&, because the internal logic of ei_nested determined that since a was already a matrix, there was no point
* in copying it into another matrix.
*/
-template<typename T, int n=1, typename PlainMatrixType = typename ei_eval<T>::type> struct ei_nested
+template<typename T, int n=1, typename PlainObject = typename ei_eval<T>::type> struct ei_nested
{
enum {
CostEval = (n+1) * int(NumTraits<typename ei_traits<T>::Scalar>::ReadCost),
@@ -266,7 +266,7 @@ template<typename T, int n=1, typename PlainMatrixType = typename ei_eval<T>::ty
typedef typename ei_meta_if<
( int(ei_traits<T>::Flags) & EvalBeforeNestingBit ) ||
( int(CostEval) <= int(CostNoEval) ),
- PlainMatrixType,
+ PlainObject,
typename ei_ref_selector<T>::type
>::ret type;
};