aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
authorGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2008-04-10 09:01:28 +0000
committerGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2008-04-10 09:01:28 +0000
commit9d8876ce82c51f7898ddf73984fcfbb413ef851f (patch)
tree6496601bb26621eb33bb37150a93db694bd9e644
parent212da8ffe0d3bad19479ffa2c5b36ce4edcd97da (diff)
* rename XprCopy -> Nested
* rename OperatorEquals -> Assign * move Util.h and FwDecl.h to a util/ subdir
-rw-r--r--Eigen/Core15
-rw-r--r--Eigen/src/Core/Assign.h (renamed from Eigen/src/Core/OperatorEquals.h)6
-rw-r--r--Eigen/src/Core/Block.h2
-rw-r--r--Eigen/src/Core/CMakeLists.txt2
-rw-r--r--Eigen/src/Core/CwiseBinaryOp.h4
-rw-r--r--Eigen/src/Core/CwiseUnaryOp.h2
-rw-r--r--Eigen/src/Core/DiagonalCoeffs.h2
-rw-r--r--Eigen/src/Core/DiagonalMatrix.h2
-rw-r--r--Eigen/src/Core/Dot.h40
-rw-r--r--Eigen/src/Core/Fuzzy.h18
-rw-r--r--Eigen/src/Core/Lazy.h2
-rw-r--r--Eigen/src/Core/Minor.h2
-rw-r--r--Eigen/src/Core/Product.h30
-rw-r--r--Eigen/src/Core/Redux.h16
-rw-r--r--Eigen/src/Core/Transpose.h2
-rw-r--r--Eigen/src/Core/util/CMakeLists.txt6
-rw-r--r--Eigen/src/Core/util/ForwardDeclarations.h (renamed from Eigen/src/Core/ForwardDeclarations.h)2
-rw-r--r--Eigen/src/Core/util/Util.h (renamed from Eigen/src/Core/Util.h)11
-rw-r--r--bench/benchmark.cpp8
19 files changed, 92 insertions, 80 deletions
diff --git a/Eigen/Core b/Eigen/Core
index dc28951f0..404eb89de 100644
--- a/Eigen/Core
+++ b/Eigen/Core
@@ -7,26 +7,17 @@
#include <cassert>
#include <iostream>
-#ifdef EIGEN_VECTORIZE
-#ifdef EIGEN_INTEL_PLATFORM
-#include <emmintrin.h>
-#include <xmmintrin.h>
-#else
-#undef EIGEN_VECTORIZE
-#endif
-#endif
-
namespace Eigen {
-#include "src/Core/Util.h"
-#include "src/Core/ForwardDeclarations.h"
+#include "src/Core/util/Util.h"
+#include "src/Core/util/ForwardDeclarations.h"
#include "src/Core/NumTraits.h"
#include "src/Core/MathFunctions.h"
#include "src/Core/PacketMath.h"
#include "src/Core/Functors.h"
#include "src/Core/MatrixBase.h"
#include "src/Core/Coeffs.h"
-#include "src/Core/OperatorEquals.h"
+#include "src/Core/Assign.h"
#include "src/Core/MatrixStorage.h"
#include "src/Core/Matrix.h"
#include "src/Core/Lazy.h"
diff --git a/Eigen/src/Core/OperatorEquals.h b/Eigen/src/Core/Assign.h
index 213c5f660..f33f491bf 100644
--- a/Eigen/src/Core/OperatorEquals.h
+++ b/Eigen/src/Core/Assign.h
@@ -23,8 +23,8 @@
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.
-#ifndef EIGEN_OPERATOREQUALS_H
-#define EIGEN_OPERATOREQUALS_H
+#ifndef EIGEN_ASSIGN_H
+#define EIGEN_ASSIGN_H
template<typename Derived1, typename Derived2, int UnrollCount>
struct ei_matrix_operator_equals_unroller
@@ -243,4 +243,4 @@ struct ei_operator_equals_impl<Derived, OtherDerived, true>
}
};
-#endif // EIGEN_OPERATOREQUALS_H
+#endif // EIGEN_ASSIGN_H
diff --git a/Eigen/src/Core/Block.h b/Eigen/src/Core/Block.h
index f0c1d11c0..c8f4ac6b1 100644
--- a/Eigen/src/Core/Block.h
+++ b/Eigen/src/Core/Block.h
@@ -143,7 +143,7 @@ template<typename MatrixType, int BlockRows, int BlockCols> class Block
protected:
- const typename MatrixType::XprCopy m_matrix;
+ const typename MatrixType::Nested m_matrix;
ei_int_if_dynamic<MatrixType::RowsAtCompileTime == 1 ? 0 : Dynamic> m_startRow;
ei_int_if_dynamic<MatrixType::ColsAtCompileTime == 1 ? 0 : Dynamic> m_startCol;
ei_int_if_dynamic<RowsAtCompileTime> m_blockRows;
diff --git a/Eigen/src/Core/CMakeLists.txt b/Eigen/src/Core/CMakeLists.txt
index 963cd50e6..a59a0fd67 100644
--- a/Eigen/src/Core/CMakeLists.txt
+++ b/Eigen/src/Core/CMakeLists.txt
@@ -4,3 +4,5 @@ INSTALL(FILES
${Eigen_Core_SRCS}
DESTINATION ${INCLUDE_INSTALL_DIR}/Eigen/src/Core
)
+
+ADD_SUBDIRECTORY(util)
diff --git a/Eigen/src/Core/CwiseBinaryOp.h b/Eigen/src/Core/CwiseBinaryOp.h
index 174b448fc..0cddc5a30 100644
--- a/Eigen/src/Core/CwiseBinaryOp.h
+++ b/Eigen/src/Core/CwiseBinaryOp.h
@@ -97,8 +97,8 @@ class CwiseBinaryOp : ei_no_assignment_operator,
}
protected:
- const typename Lhs::XprCopy m_lhs;
- const typename Rhs::XprCopy m_rhs;
+ const typename Lhs::Nested m_lhs;
+ const typename Rhs::Nested m_rhs;
const BinaryOp m_functor;
};
diff --git a/Eigen/src/Core/CwiseUnaryOp.h b/Eigen/src/Core/CwiseUnaryOp.h
index cbe545a92..2a6858703 100644
--- a/Eigen/src/Core/CwiseUnaryOp.h
+++ b/Eigen/src/Core/CwiseUnaryOp.h
@@ -83,7 +83,7 @@ class CwiseUnaryOp : ei_no_assignment_operator,
}
protected:
- const typename MatrixType::XprCopy m_matrix;
+ const typename MatrixType::Nested m_matrix;
const UnaryOp m_functor;
};
diff --git a/Eigen/src/Core/DiagonalCoeffs.h b/Eigen/src/Core/DiagonalCoeffs.h
index 7f8fea162..3e4ebded8 100644
--- a/Eigen/src/Core/DiagonalCoeffs.h
+++ b/Eigen/src/Core/DiagonalCoeffs.h
@@ -87,7 +87,7 @@ template<typename MatrixType> class DiagonalCoeffs
protected:
- const typename MatrixType::XprCopy m_matrix;
+ const typename MatrixType::Nested m_matrix;
};
/** \returns an expression of the main diagonal of the matrix \c *this
diff --git a/Eigen/src/Core/DiagonalMatrix.h b/Eigen/src/Core/DiagonalMatrix.h
index 6a243a402..0fae56805 100644
--- a/Eigen/src/Core/DiagonalMatrix.h
+++ b/Eigen/src/Core/DiagonalMatrix.h
@@ -77,7 +77,7 @@ class DiagonalMatrix : ei_no_assignment_operator,
}
protected:
- const typename CoeffsVectorType::XprCopy m_coeffs;
+ const typename CoeffsVectorType::Nested m_coeffs;
};
/** \returns an expression of a diagonal matrix with *this as vector of diagonal coefficients
diff --git a/Eigen/src/Core/Dot.h b/Eigen/src/Core/Dot.h
index 1d768b259..ae42bfa1b 100644
--- a/Eigen/src/Core/Dot.h
+++ b/Eigen/src/Core/Dot.h
@@ -72,31 +72,31 @@ template<typename OtherDerived>
typename ei_traits<Derived>::Scalar
MatrixBase<Derived>::dot(const MatrixBase<OtherDerived>& other) const
{
- typedef typename Derived::XprCopy XprCopy;
- typedef typename OtherDerived::XprCopy OtherXprCopy;
- typedef typename ei_unref<XprCopy>::type _XprCopy;
- typedef typename ei_unref<OtherXprCopy>::type _OtherXprCopy;
- XprCopy xprCopy(derived());
- OtherXprCopy otherXprCopy(other.derived());
+ typedef typename Derived::Nested Nested;
+ typedef typename OtherDerived::Nested OtherNested;
+ typedef typename ei_unref<Nested>::type _Nested;
+ typedef typename ei_unref<OtherNested>::type _OtherNested;
+ Nested nested(derived());
+ OtherNested otherNested(other.derived());
- ei_assert(_XprCopy::IsVectorAtCompileTime
- && _OtherXprCopy::IsVectorAtCompileTime
- && xprCopy.size() == otherXprCopy.size());
+ ei_assert(_Nested::IsVectorAtCompileTime
+ && _OtherNested::IsVectorAtCompileTime
+ && nested.size() == otherNested.size());
Scalar res;
const bool unroll = SizeAtCompileTime
- * (_XprCopy::CoeffReadCost + _OtherXprCopy::CoeffReadCost + NumTraits<Scalar>::MulCost)
+ * (_Nested::CoeffReadCost + _OtherNested::CoeffReadCost + NumTraits<Scalar>::MulCost)
+ (SizeAtCompileTime - 1) * NumTraits<Scalar>::AddCost
<= EIGEN_UNROLLING_LIMIT;
if(unroll)
ei_dot_unroller<SizeAtCompileTime-1,
unroll ? SizeAtCompileTime : Dynamic,
- _XprCopy, _OtherXprCopy>
- ::run(xprCopy, otherXprCopy, res);
+ _Nested, _OtherNested>
+ ::run(nested, otherNested, res);
else
{
- res = xprCopy.coeff(0) * ei_conj(otherXprCopy.coeff(0));
+ res = nested.coeff(0) * ei_conj(otherNested.coeff(0));
for(int i = 1; i < size(); i++)
- res += xprCopy.coeff(i)* ei_conj(otherXprCopy.coeff(i));
+ res += nested.coeff(i)* ei_conj(otherNested.coeff(i));
}
return res;
}
@@ -149,9 +149,9 @@ template<typename OtherDerived>
bool MatrixBase<Derived>::isOrtho
(const MatrixBase<OtherDerived>& other, RealScalar prec) const
{
- typename ei_xpr_copy<Derived,2>::type xprCopy(derived());
- typename ei_xpr_copy<OtherDerived,2>::type otherXprCopy(other.derived());
- return ei_abs2(xprCopy.dot(otherXprCopy)) <= prec * prec * xprCopy.norm2() * otherXprCopy.norm2();
+ typename ei_nested<Derived,2>::type nested(derived());
+ typename ei_nested<OtherDerived,2>::type otherNested(other.derived());
+ return ei_abs2(nested.dot(otherNested)) <= prec * prec * nested.norm2() * otherNested.norm2();
}
/** \returns true if *this is approximately an unitary matrix,
@@ -168,13 +168,13 @@ bool MatrixBase<Derived>::isOrtho
template<typename Derived>
bool MatrixBase<Derived>::isOrtho(RealScalar prec) const
{
- typename Derived::XprCopy xprCopy(derived());
+ typename Derived::Nested nested(derived());
for(int i = 0; i < cols(); i++)
{
- if(!ei_isApprox(xprCopy.col(i).norm2(), static_cast<Scalar>(1), prec))
+ if(!ei_isApprox(nested.col(i).norm2(), static_cast<Scalar>(1), prec))
return false;
for(int j = 0; j < i; j++)
- if(!ei_isMuchSmallerThan(xprCopy.col(i).dot(xprCopy.col(j)), static_cast<Scalar>(1), prec))
+ if(!ei_isMuchSmallerThan(nested.col(i).dot(nested.col(j)), static_cast<Scalar>(1), prec))
return false;
}
return true;
diff --git a/Eigen/src/Core/Fuzzy.h b/Eigen/src/Core/Fuzzy.h
index 5dd0265ba..8438193d8 100644
--- a/Eigen/src/Core/Fuzzy.h
+++ b/Eigen/src/Core/Fuzzy.h
@@ -55,11 +55,11 @@ bool MatrixBase<Derived>::isApprox(
}
else
{
- typename Derived::XprCopy xprCopy(derived());
- typename OtherDerived::XprCopy otherXprCopy(other.derived());
+ typename Derived::Nested nested(derived());
+ typename OtherDerived::Nested otherNested(other.derived());
for(int i = 0; i < cols(); i++)
- if((xprCopy.col(i) - otherXprCopy.col(i)).norm2()
- > std::min(xprCopy.col(i).norm2(), otherXprCopy.col(i).norm2()) * prec * prec)
+ if((nested.col(i) - otherNested.col(i)).norm2()
+ > std::min(nested.col(i).norm2(), otherNested.col(i).norm2()) * prec * prec)
return false;
return true;
}
@@ -87,9 +87,9 @@ bool MatrixBase<Derived>::isMuchSmallerThan(
}
else
{
- typename Derived::XprCopy xprCopy(*this);
+ typename Derived::Nested nested(*this);
for(int i = 0; i < cols(); i++)
- if(xprCopy.col(i).norm2() > ei_abs2(other * prec))
+ if(nested.col(i).norm2() > ei_abs2(other * prec))
return false;
return true;
}
@@ -119,10 +119,10 @@ bool MatrixBase<Derived>::isMuchSmallerThan(
}
else
{
- typename Derived::XprCopy xprCopy(*this);
- typename OtherDerived::XprCopy otherXprCopy(other);
+ typename Derived::Nested nested(*this);
+ typename OtherDerived::Nested otherNested(other);
for(int i = 0; i < cols(); i++)
- if(xprCopy.col(i).norm2() > otherXprCopy.col(i).norm2() * prec * prec)
+ if(nested.col(i).norm2() > otherNested.col(i).norm2() * prec * prec)
return false;
return true;
}
diff --git a/Eigen/src/Core/Lazy.h b/Eigen/src/Core/Lazy.h
index aacc61695..0c65cdeba 100644
--- a/Eigen/src/Core/Lazy.h
+++ b/Eigen/src/Core/Lazy.h
@@ -73,7 +73,7 @@ template<typename ExpressionType> class Lazy
}
protected:
- const typename ExpressionType::XprCopy m_expression;
+ const typename ExpressionType::Nested m_expression;
};
/** \returns an expression of the lazy version of *this.
diff --git a/Eigen/src/Core/Minor.h b/Eigen/src/Core/Minor.h
index 1b060928f..5c9afc162 100644
--- a/Eigen/src/Core/Minor.h
+++ b/Eigen/src/Core/Minor.h
@@ -88,7 +88,7 @@ template<typename MatrixType> class Minor
}
protected:
- const typename MatrixType::XprCopy m_matrix;
+ const typename MatrixType::Nested m_matrix;
const int m_row, m_col;
};
diff --git a/Eigen/src/Core/Product.h b/Eigen/src/Core/Product.h
index 142955108..9e7127720 100644
--- a/Eigen/src/Core/Product.h
+++ b/Eigen/src/Core/Product.h
@@ -116,15 +116,15 @@ template<typename Lhs, typename Rhs, int EvalMode>
struct ei_traits<Product<Lhs, Rhs, EvalMode> >
{
typedef typename Lhs::Scalar Scalar;
- typedef typename ei_xpr_copy<Lhs,Rhs::ColsAtCompileTime>::type LhsXprCopy;
- typedef typename ei_xpr_copy<Rhs,Lhs::RowsAtCompileTime>::type RhsXprCopy;
- typedef typename ei_unref<LhsXprCopy>::type _LhsXprCopy;
- typedef typename ei_unref<RhsXprCopy>::type _RhsXprCopy;
+ typedef typename ei_nested<Lhs,Rhs::ColsAtCompileTime>::type LhsNested;
+ typedef typename ei_nested<Rhs,Lhs::RowsAtCompileTime>::type RhsNested;
+ typedef typename ei_unref<LhsNested>::type _LhsNested;
+ typedef typename ei_unref<RhsNested>::type _RhsNested;
enum {
- LhsCoeffReadCost = _LhsXprCopy::CoeffReadCost,
- RhsCoeffReadCost = _RhsXprCopy::CoeffReadCost,
- LhsFlags = _LhsXprCopy::Flags,
- RhsFlags = _RhsXprCopy::Flags,
+ LhsCoeffReadCost = _LhsNested::CoeffReadCost,
+ RhsCoeffReadCost = _RhsNested::CoeffReadCost,
+ LhsFlags = _LhsNested::Flags,
+ RhsFlags = _RhsNested::Flags,
RowsAtCompileTime = Lhs::RowsAtCompileTime,
ColsAtCompileTime = Rhs::ColsAtCompileTime,
MaxRowsAtCompileTime = Lhs::MaxRowsAtCompileTime,
@@ -153,10 +153,10 @@ template<typename Lhs, typename Rhs, int EvalMode> class Product : ei_no_assignm
public:
EIGEN_GENERIC_PUBLIC_INTERFACE(Product)
- typedef typename ei_traits<Product>::LhsXprCopy LhsXprCopy;
- typedef typename ei_traits<Product>::RhsXprCopy RhsXprCopy;
- typedef typename ei_traits<Product>::_LhsXprCopy _LhsXprCopy;
- typedef typename ei_traits<Product>::_RhsXprCopy _RhsXprCopy;
+ typedef typename ei_traits<Product>::LhsNested LhsNested;
+ typedef typename ei_traits<Product>::RhsNested RhsNested;
+ typedef typename ei_traits<Product>::_LhsNested _LhsNested;
+ typedef typename ei_traits<Product>::_RhsNested _RhsNested;
Product(const Lhs& lhs, const Rhs& rhs)
: m_lhs(lhs), m_rhs(rhs)
@@ -181,7 +181,7 @@ template<typename Lhs, typename Rhs, int EvalMode> class Product : ei_no_assignm
{
ei_product_unroller<Lhs::ColsAtCompileTime-1,
unroll ? Lhs::ColsAtCompileTime : Dynamic,
- _LhsXprCopy, _RhsXprCopy>
+ _LhsNested, _RhsNested>
::run(row, col, m_lhs, m_rhs, res);
}
else
@@ -224,8 +224,8 @@ template<typename Lhs, typename Rhs, int EvalMode> class Product : ei_no_assignm
}
protected:
- const LhsXprCopy m_lhs;
- const RhsXprCopy m_rhs;
+ const LhsNested m_lhs;
+ const RhsNested m_rhs;
};
/** \returns the matrix product of \c *this and \a other.
diff --git a/Eigen/src/Core/Redux.h b/Eigen/src/Core/Redux.h
index 12ceedd76..6e740fc7d 100644
--- a/Eigen/src/Core/Redux.h
+++ b/Eigen/src/Core/Redux.h
@@ -87,18 +87,18 @@ struct ei_traits<PartialRedux<Direction, BinaryOp, MatrixType> >
typedef typename ei_result_of<
BinaryOp(typename MatrixType::Scalar)
>::type Scalar;
- typedef typename ei_xpr_copy<MatrixType>::type MatrixTypeXprCopy;
- typedef typename ei_unref<MatrixTypeXprCopy>::type _MatrixTypeXprCopy;
+ typedef typename ei_nested<MatrixType>::type MatrixTypeNested;
+ typedef typename ei_unref<MatrixTypeNested>::type _MatrixTypeNested;
enum {
RowsAtCompileTime = Direction==Vertical ? 1 : MatrixType::RowsAtCompileTime,
ColsAtCompileTime = Direction==Horizontal ? 1 : MatrixType::ColsAtCompileTime,
MaxRowsAtCompileTime = Direction==Vertical ? 1 : MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = Direction==Horizontal ? 1 : MatrixType::MaxColsAtCompileTime,
Flags = ((RowsAtCompileTime == Dynamic || ColsAtCompileTime == Dynamic)
- ? (unsigned int)_MatrixTypeXprCopy::Flags
- : (unsigned int)_MatrixTypeXprCopy::Flags & ~LargeBit) & ~VectorizableBit,
+ ? (unsigned int)_MatrixTypeNested::Flags
+ : (unsigned int)_MatrixTypeNested::Flags & ~LargeBit) & ~VectorizableBit,
TraversalSize = Direction==Vertical ? RowsAtCompileTime : ColsAtCompileTime,
- CoeffReadCost = TraversalSize * _MatrixTypeXprCopy::CoeffReadCost
+ CoeffReadCost = TraversalSize * _MatrixTypeNested::CoeffReadCost
+ (TraversalSize - 1) * ei_functor_traits<BinaryOp>::Cost
};
};
@@ -110,8 +110,8 @@ class PartialRedux : ei_no_assignment_operator,
public:
EIGEN_GENERIC_PUBLIC_INTERFACE(PartialRedux)
- typedef typename ei_traits<PartialRedux>::MatrixTypeXprCopy MatrixTypeXprCopy;
- typedef typename ei_traits<PartialRedux>::_MatrixTypeXprCopy _MatrixTypeXprCopy;
+ typedef typename ei_traits<PartialRedux>::MatrixTypeNested MatrixTypeNested;
+ typedef typename ei_traits<PartialRedux>::_MatrixTypeNested _MatrixTypeNested;
PartialRedux(const MatrixType& mat, const BinaryOp& func = BinaryOp())
: m_matrix(mat), m_functor(func) {}
@@ -130,7 +130,7 @@ class PartialRedux : ei_no_assignment_operator,
}
protected:
- const MatrixTypeXprCopy m_matrix;
+ const MatrixTypeNested m_matrix;
const BinaryOp m_functor;
};
diff --git a/Eigen/src/Core/Transpose.h b/Eigen/src/Core/Transpose.h
index 6710f3092..424557754 100644
--- a/Eigen/src/Core/Transpose.h
+++ b/Eigen/src/Core/Transpose.h
@@ -88,7 +88,7 @@ template<typename MatrixType> class Transpose
}
protected:
- const typename MatrixType::XprCopy m_matrix;
+ const typename MatrixType::Nested m_matrix;
};
/** \returns an expression of the transpose of *this.
diff --git a/Eigen/src/Core/util/CMakeLists.txt b/Eigen/src/Core/util/CMakeLists.txt
new file mode 100644
index 000000000..9ab9fc8b0
--- /dev/null
+++ b/Eigen/src/Core/util/CMakeLists.txt
@@ -0,0 +1,6 @@
+FILE(GLOB Eigen_Core_util_SRCS "*.h")
+
+INSTALL(FILES
+ ${Eigen_Core_SRCS}
+ DESTINATION ${INCLUDE_INSTALL_DIR}/Eigen/src/Core/util
+ )
diff --git a/Eigen/src/Core/ForwardDeclarations.h b/Eigen/src/Core/util/ForwardDeclarations.h
index d9699301d..a7f269cb4 100644
--- a/Eigen/src/Core/ForwardDeclarations.h
+++ b/Eigen/src/Core/util/ForwardDeclarations.h
@@ -93,7 +93,7 @@ template<typename T> struct ei_is_temporary<Temporary<T> >
enum { ret = 1 };
};
-template<typename T, int n=1> struct ei_xpr_copy
+template<typename T, int n=1> struct ei_nested
{
typedef typename ei_meta_if<
ei_is_temporary<T>::ret,
diff --git a/Eigen/src/Core/Util.h b/Eigen/src/Core/util/Util.h
index e2c95bc53..de45927ce 100644
--- a/Eigen/src/Core/Util.h
+++ b/Eigen/src/Core/util/Util.h
@@ -26,6 +26,15 @@
#ifndef EIGEN_UTIL_H
#define EIGEN_UTIL_H
+#ifdef EIGEN_VECTORIZE
+#ifdef EIGEN_INTEL_PLATFORM
+#include <emmintrin.h>
+#include <xmmintrin.h>
+#else
+#undef EIGEN_VECTORIZE
+#endif
+#endif
+
#ifdef EIGEN_DONT_USE_UNROLLED_LOOPS
#define EIGEN_UNROLLING_LIMIT 0
#endif
@@ -117,7 +126,7 @@ EIGEN_INHERIT_SCALAR_ASSIGNMENT_OPERATOR(Derived, /=)
typedef BaseClass Base; \
typedef typename Eigen::ei_traits<Derived>::Scalar Scalar; \
typedef typename Base::PacketScalar PacketScalar; \
-typedef typename Eigen::ei_xpr_copy<Derived>::type XprCopy; \
+typedef typename Eigen::ei_nested<Derived>::type Nested; \
typedef typename Eigen::ei_eval<Derived>::type Eval; \
enum { RowsAtCompileTime = Base::RowsAtCompileTime, \
ColsAtCompileTime = Base::ColsAtCompileTime, \
diff --git a/bench/benchmark.cpp b/bench/benchmark.cpp
index 4ff678d8a..99fd097ba 100644
--- a/bench/benchmark.cpp
+++ b/bench/benchmark.cpp
@@ -12,10 +12,14 @@ USING_PART_OF_NAMESPACE_EIGEN
#define REPEAT 40000000
#endif
+#ifndef SCALAR
+#define SCALAR double
+#endif
+
int main(int argc, char *argv[])
{
- Matrix<double,MATSIZE,MATSIZE> I;
- Matrix<double,MATSIZE,MATSIZE> m;
+ Matrix<SCALAR,MATSIZE,MATSIZE> I;
+ Matrix<SCALAR,MATSIZE,MATSIZE> m;
for(int i = 0; i < MATSIZE; i++)
for(int j = 0; j < MATSIZE; j++)
{