aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2015-06-19 17:56:39 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2015-06-19 17:56:39 +0200
commit846b227bb7cb75b15d56b864fd294e028a33db3a (patch)
treee43c0703640ad6767bf8707972632a1b808ba6ed
parent386d9e5ebdff7f986a4f707e965703f2785dc292 (diff)
Get rid of class internal::nested<> (still have to updated Tensor module)
-rw-r--r--Eigen/src/Core/ArrayWrapper.h4
-rw-r--r--Eigen/src/Core/Block.h2
-rw-r--r--Eigen/src/Core/CwiseBinaryOp.h4
-rw-r--r--Eigen/src/Core/CwiseUnaryView.h3
-rw-r--r--Eigen/src/Core/Diagonal.h2
-rw-r--r--Eigen/src/Core/Inverse.h2
-rw-r--r--Eigen/src/Core/Product.h4
-rw-r--r--Eigen/src/Core/Replicate.h5
-rw-r--r--Eigen/src/Core/Reverse.h2
-rw-r--r--Eigen/src/Core/SelfAdjointView.h2
-rw-r--r--Eigen/src/Core/Transpose.h2
-rw-r--r--Eigen/src/Core/TriangularMatrix.h2
-rw-r--r--Eigen/src/Core/VectorwiseOp.h2
-rw-r--r--Eigen/src/Core/util/Macros.h2
-rw-r--r--Eigen/src/Core/util/XprHelper.h9
-rw-r--r--Eigen/src/Geometry/Homogeneous.h2
-rw-r--r--Eigen/src/SparseCore/SparseMatrix.h2
-rw-r--r--Eigen/src/SparseCore/SparseUtil.h2
-rw-r--r--unsupported/Eigen/CXX11/src/Tensor/TensorTraits.h5
-rw-r--r--unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h2
-rw-r--r--unsupported/Eigen/src/MatrixFunctions/MatrixFunction.h2
-rw-r--r--unsupported/Eigen/src/MatrixFunctions/MatrixLogarithm.h2
-rw-r--r--unsupported/Eigen/src/MatrixFunctions/MatrixSquareRoot.h2
-rw-r--r--unsupported/Eigen/src/Skyline/SkylineProduct.h4
-rw-r--r--unsupported/Eigen/src/SparseExtra/BlockSparseMatrix.h2
25 files changed, 33 insertions, 39 deletions
diff --git a/Eigen/src/Core/ArrayWrapper.h b/Eigen/src/Core/ArrayWrapper.h
index 2ef112986..4e484f290 100644
--- a/Eigen/src/Core/ArrayWrapper.h
+++ b/Eigen/src/Core/ArrayWrapper.h
@@ -52,7 +52,7 @@ class ArrayWrapper : public ArrayBase<ArrayWrapper<ExpressionType> >
const Scalar
>::type ScalarWithConstIfNotLvalue;
- typedef typename internal::nested<ExpressionType>::type NestedExpressionType;
+ typedef typename internal::ref_selector<ExpressionType>::type NestedExpressionType;
EIGEN_DEVICE_FUNC
explicit EIGEN_STRONG_INLINE ArrayWrapper(ExpressionType& matrix) : m_expression(matrix) {}
@@ -195,7 +195,7 @@ class MatrixWrapper : public MatrixBase<MatrixWrapper<ExpressionType> >
const Scalar
>::type ScalarWithConstIfNotLvalue;
- typedef typename internal::nested<ExpressionType>::type NestedExpressionType;
+ typedef typename internal::ref_selector<ExpressionType>::type NestedExpressionType;
EIGEN_DEVICE_FUNC
explicit inline MatrixWrapper(ExpressionType& matrix) : m_expression(matrix) {}
diff --git a/Eigen/src/Core/Block.h b/Eigen/src/Core/Block.h
index 7f84534e1..aed6147c7 100644
--- a/Eigen/src/Core/Block.h
+++ b/Eigen/src/Core/Block.h
@@ -55,7 +55,7 @@ struct traits<Block<XprType, BlockRows, BlockCols, InnerPanel> > : traits<XprTyp
typedef typename traits<XprType>::Scalar Scalar;
typedef typename traits<XprType>::StorageKind StorageKind;
typedef typename traits<XprType>::XprKind XprKind;
- typedef typename nested<XprType>::type XprTypeNested;
+ typedef typename ref_selector<XprType>::type XprTypeNested;
typedef typename remove_reference<XprTypeNested>::type _XprTypeNested;
enum{
MatrixRows = traits<XprType>::RowsAtCompileTime,
diff --git a/Eigen/src/Core/CwiseBinaryOp.h b/Eigen/src/Core/CwiseBinaryOp.h
index 0c9d9fbf2..e42c3031b 100644
--- a/Eigen/src/Core/CwiseBinaryOp.h
+++ b/Eigen/src/Core/CwiseBinaryOp.h
@@ -95,8 +95,8 @@ class CwiseBinaryOp :
BinaryOp>::ret>::Base Base;
EIGEN_GENERIC_PUBLIC_INTERFACE(CwiseBinaryOp)
- typedef typename internal::nested<LhsType>::type LhsNested;
- typedef typename internal::nested<RhsType>::type RhsNested;
+ typedef typename internal::ref_selector<LhsType>::type LhsNested;
+ typedef typename internal::ref_selector<RhsType>::type RhsNested;
typedef typename internal::remove_reference<LhsNested>::type _LhsNested;
typedef typename internal::remove_reference<RhsNested>::type _RhsNested;
diff --git a/Eigen/src/Core/CwiseUnaryView.h b/Eigen/src/Core/CwiseUnaryView.h
index 6680f32dd..72244751e 100644
--- a/Eigen/src/Core/CwiseUnaryView.h
+++ b/Eigen/src/Core/CwiseUnaryView.h
@@ -84,8 +84,7 @@ class CwiseUnaryView : public CwiseUnaryViewImpl<ViewOp, MatrixType, typename in
nestedExpression() { return m_matrix.const_cast_derived(); }
protected:
- // FIXME changed from MatrixType::Nested because of a weird compilation error with sun CC
- typename internal::nested<MatrixType>::type m_matrix;
+ typename internal::ref_selector<MatrixType>::type m_matrix;
ViewOp m_functor;
};
diff --git a/Eigen/src/Core/Diagonal.h b/Eigen/src/Core/Diagonal.h
index 2446a18d4..c79891ae1 100644
--- a/Eigen/src/Core/Diagonal.h
+++ b/Eigen/src/Core/Diagonal.h
@@ -37,7 +37,7 @@ template<typename MatrixType, int DiagIndex>
struct traits<Diagonal<MatrixType,DiagIndex> >
: traits<MatrixType>
{
- typedef typename nested<MatrixType>::type MatrixTypeNested;
+ typedef typename ref_selector<MatrixType>::type MatrixTypeNested;
typedef typename remove_reference<MatrixTypeNested>::type _MatrixTypeNested;
typedef typename MatrixType::StorageKind StorageKind;
enum {
diff --git a/Eigen/src/Core/Inverse.h b/Eigen/src/Core/Inverse.h
index f3fa82a01..53d9b28f5 100644
--- a/Eigen/src/Core/Inverse.h
+++ b/Eigen/src/Core/Inverse.h
@@ -47,7 +47,7 @@ class Inverse : public InverseImpl<XprType,typename internal::traits<XprType>::S
public:
typedef typename XprType::StorageIndex StorageIndex;
typedef typename XprType::PlainObject PlainObject;
- typedef typename internal::nested<XprType>::type XprTypeNested;
+ typedef typename internal::ref_selector<XprType>::type XprTypeNested;
typedef typename internal::remove_all<XprTypeNested>::type XprTypeNestedCleaned;
explicit Inverse(const XprType &xpr)
diff --git a/Eigen/src/Core/Product.h b/Eigen/src/Core/Product.h
index 87340b88c..f58fa5f60 100644
--- a/Eigen/src/Core/Product.h
+++ b/Eigen/src/Core/Product.h
@@ -121,8 +121,8 @@ class Product : public ProductImpl<_Lhs,_Rhs,Option,
internal::product_type<Lhs,Rhs>::ret>::ret>::Base Base;
EIGEN_GENERIC_PUBLIC_INTERFACE(Product)
- typedef typename internal::nested<Lhs>::type LhsNested;
- typedef typename internal::nested<Rhs>::type RhsNested;
+ typedef typename internal::ref_selector<Lhs>::type LhsNested;
+ typedef typename internal::ref_selector<Rhs>::type RhsNested;
typedef typename internal::remove_all<LhsNested>::type LhsNestedCleaned;
typedef typename internal::remove_all<RhsNested>::type RhsNestedCleaned;
diff --git a/Eigen/src/Core/Replicate.h b/Eigen/src/Core/Replicate.h
index 208f86380..d105de62b 100644
--- a/Eigen/src/Core/Replicate.h
+++ b/Eigen/src/Core/Replicate.h
@@ -35,10 +35,7 @@ struct traits<Replicate<MatrixType,RowFactor,ColFactor> >
typedef typename MatrixType::Scalar Scalar;
typedef typename traits<MatrixType>::StorageKind StorageKind;
typedef typename traits<MatrixType>::XprKind XprKind;
- enum {
- Factor = (RowFactor==Dynamic || ColFactor==Dynamic) ? Dynamic : RowFactor*ColFactor
- };
- typedef typename nested<MatrixType,Factor>::type MatrixTypeNested;
+ typedef typename ref_selector<MatrixType>::type MatrixTypeNested;
typedef typename remove_reference<MatrixTypeNested>::type _MatrixTypeNested;
enum {
RowsAtCompileTime = RowFactor==Dynamic || int(MatrixType::RowsAtCompileTime)==Dynamic
diff --git a/Eigen/src/Core/Reverse.h b/Eigen/src/Core/Reverse.h
index 5237fbf1c..993458aa5 100644
--- a/Eigen/src/Core/Reverse.h
+++ b/Eigen/src/Core/Reverse.h
@@ -37,7 +37,7 @@ struct traits<Reverse<MatrixType, Direction> >
typedef typename MatrixType::Scalar Scalar;
typedef typename traits<MatrixType>::StorageKind StorageKind;
typedef typename traits<MatrixType>::XprKind XprKind;
- typedef typename nested<MatrixType>::type MatrixTypeNested;
+ typedef typename ref_selector<MatrixType>::type MatrixTypeNested;
typedef typename remove_reference<MatrixTypeNested>::type _MatrixTypeNested;
enum {
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
diff --git a/Eigen/src/Core/SelfAdjointView.h b/Eigen/src/Core/SelfAdjointView.h
index a05746ad2..87e87ab3a 100644
--- a/Eigen/src/Core/SelfAdjointView.h
+++ b/Eigen/src/Core/SelfAdjointView.h
@@ -32,7 +32,7 @@ namespace internal {
template<typename MatrixType, unsigned int UpLo>
struct traits<SelfAdjointView<MatrixType, UpLo> > : traits<MatrixType>
{
- typedef typename nested<MatrixType>::type MatrixTypeNested;
+ typedef typename ref_selector<MatrixType>::type MatrixTypeNested;
typedef typename remove_all<MatrixTypeNested>::type MatrixTypeNestedCleaned;
typedef MatrixType ExpressionType;
typedef typename MatrixType::PlainObject FullMatrixType;
diff --git a/Eigen/src/Core/Transpose.h b/Eigen/src/Core/Transpose.h
index 576411600..e205cec4a 100644
--- a/Eigen/src/Core/Transpose.h
+++ b/Eigen/src/Core/Transpose.h
@@ -31,7 +31,7 @@ namespace internal {
template<typename MatrixType>
struct traits<Transpose<MatrixType> > : public traits<MatrixType>
{
- typedef typename nested<MatrixType>::type MatrixTypeNested;
+ typedef typename ref_selector<MatrixType>::type MatrixTypeNested;
typedef typename remove_reference<MatrixTypeNested>::type MatrixTypeNestedPlain;
enum {
RowsAtCompileTime = MatrixType::ColsAtCompileTime,
diff --git a/Eigen/src/Core/TriangularMatrix.h b/Eigen/src/Core/TriangularMatrix.h
index 0a6397509..bf03ea818 100644
--- a/Eigen/src/Core/TriangularMatrix.h
+++ b/Eigen/src/Core/TriangularMatrix.h
@@ -164,7 +164,7 @@ namespace internal {
template<typename MatrixType, unsigned int _Mode>
struct traits<TriangularView<MatrixType, _Mode> > : traits<MatrixType>
{
- typedef typename nested<MatrixType>::type MatrixTypeNested;
+ typedef typename ref_selector<MatrixType>::type MatrixTypeNested;
typedef typename remove_reference<MatrixTypeNested>::type MatrixTypeNestedNonRef;
typedef typename remove_all<MatrixTypeNested>::type MatrixTypeNestedCleaned;
typedef typename MatrixType::PlainObject FullMatrixType;
diff --git a/Eigen/src/Core/VectorwiseOp.h b/Eigen/src/Core/VectorwiseOp.h
index ea3d8f4b1..a0f10b17b 100644
--- a/Eigen/src/Core/VectorwiseOp.h
+++ b/Eigen/src/Core/VectorwiseOp.h
@@ -41,7 +41,7 @@ struct traits<PartialReduxExpr<MatrixType, MemberOp, Direction> >
typedef typename traits<MatrixType>::StorageKind StorageKind;
typedef typename traits<MatrixType>::XprKind XprKind;
typedef typename MatrixType::Scalar InputScalar;
- typedef typename nested<MatrixType>::type MatrixTypeNested;
+ typedef typename ref_selector<MatrixType>::type MatrixTypeNested;
typedef typename remove_all<MatrixTypeNested>::type _MatrixTypeNested;
enum {
RowsAtCompileTime = Direction==Vertical ? 1 : MatrixType::RowsAtCompileTime,
diff --git a/Eigen/src/Core/util/Macros.h b/Eigen/src/Core/util/Macros.h
index 583a0cbd2..8efe2f4ca 100644
--- a/Eigen/src/Core/util/Macros.h
+++ b/Eigen/src/Core/util/Macros.h
@@ -685,7 +685,7 @@ namespace Eigen {
typedef typename Eigen::internal::traits<Derived>::Scalar Scalar; /*!< \brief Numeric type, e.g. float, double, int or std::complex<float>. */ \
typedef typename Eigen::NumTraits<Scalar>::Real RealScalar; /*!< \brief The underlying numeric type for composed scalar types. \details In cases where Scalar is e.g. std::complex<T>, T were corresponding to RealScalar. */ \
typedef typename Base::CoeffReturnType CoeffReturnType; /*!< \brief The return type for coefficient access. \details Depending on whether the object allows direct coefficient access (e.g. for a MatrixXd), this type is either 'const Scalar&' or simply 'Scalar' for objects that do not allow direct coefficient access. */ \
- typedef typename Eigen::internal::nested<Derived>::type Nested; \
+ typedef typename Eigen::internal::ref_selector<Derived>::type Nested; \
typedef typename Eigen::internal::traits<Derived>::StorageKind StorageKind; \
typedef typename Eigen::internal::traits<Derived>::StorageIndex StorageIndex; \
enum { RowsAtCompileTime = Eigen::internal::traits<Derived>::RowsAtCompileTime, \
diff --git a/Eigen/src/Core/util/XprHelper.h b/Eigen/src/Core/util/XprHelper.h
index 55132c8cf..8e12ce5a0 100644
--- a/Eigen/src/Core/util/XprHelper.h
+++ b/Eigen/src/Core/util/XprHelper.h
@@ -312,7 +312,7 @@ template<typename T> struct plain_matrix_type_row_major
> type;
};
-// we should be able to get rid of this one too
+// TODO we should be able to get rid of this one too
template<typename T> struct must_nest_by_value { enum { ret = false }; };
/** \internal The reference selector for template expressions. The idea is that we don't
@@ -340,13 +340,6 @@ struct transfer_constness
};
-// When using evaluators, we never evaluate when assembling the expression!!
-// TODO: get rid of this nested class since it's just an alias for ref_selector.
-template<typename T, int n=1, typename PlainObject = void> struct nested
-{
- typedef typename ref_selector<T>::type type;
-};
-
// However, we still need a mechanism to detect whether an expression which is evaluated multiple time
// has to be evaluated into a temporary.
// That's the purpose of this new nested_eval helper:
diff --git a/Eigen/src/Geometry/Homogeneous.h b/Eigen/src/Geometry/Homogeneous.h
index b7f996615..9f86eeb3d 100644
--- a/Eigen/src/Geometry/Homogeneous.h
+++ b/Eigen/src/Geometry/Homogeneous.h
@@ -34,7 +34,7 @@ struct traits<Homogeneous<MatrixType,Direction> >
: traits<MatrixType>
{
typedef typename traits<MatrixType>::StorageKind StorageKind;
- typedef typename nested<MatrixType>::type MatrixTypeNested;
+ typedef typename ref_selector<MatrixType>::type MatrixTypeNested;
typedef typename remove_reference<MatrixTypeNested>::type _MatrixTypeNested;
enum {
RowsPlusOne = (MatrixType::RowsAtCompileTime != Dynamic) ?
diff --git a/Eigen/src/SparseCore/SparseMatrix.h b/Eigen/src/SparseCore/SparseMatrix.h
index d9964d0f6..fb05fb4b6 100644
--- a/Eigen/src/SparseCore/SparseMatrix.h
+++ b/Eigen/src/SparseCore/SparseMatrix.h
@@ -60,7 +60,7 @@ template<typename _Scalar, int _Options, typename _Index, int DiagIndex>
struct traits<Diagonal<SparseMatrix<_Scalar, _Options, _Index>, DiagIndex> >
{
typedef SparseMatrix<_Scalar, _Options, _Index> MatrixType;
- typedef typename nested<MatrixType>::type MatrixTypeNested;
+ typedef typename ref_selector<MatrixType>::type MatrixTypeNested;
typedef typename remove_reference<MatrixTypeNested>::type _MatrixTypeNested;
typedef _Scalar Scalar;
diff --git a/Eigen/src/SparseCore/SparseUtil.h b/Eigen/src/SparseCore/SparseUtil.h
index ac2b4f10b..215601aa6 100644
--- a/Eigen/src/SparseCore/SparseUtil.h
+++ b/Eigen/src/SparseCore/SparseUtil.h
@@ -47,7 +47,7 @@ EIGEN_SPARSE_INHERIT_SCALAR_ASSIGNMENT_OPERATOR(Derived, /=)
#define _EIGEN_SPARSE_PUBLIC_INTERFACE(Derived) \
typedef typename Eigen::internal::traits<Derived >::Scalar Scalar; \
typedef typename Eigen::NumTraits<Scalar>::Real RealScalar; \
- typedef typename Eigen::internal::nested<Derived >::type Nested; \
+ typedef typename Eigen::internal::ref_selector<Derived >::type Nested; \
typedef typename Eigen::internal::traits<Derived >::StorageKind StorageKind; \
typedef typename Eigen::internal::traits<Derived >::StorageIndex StorageIndex; \
enum { RowsAtCompileTime = Eigen::internal::traits<Derived >::RowsAtCompileTime, \
diff --git a/unsupported/Eigen/CXX11/src/Tensor/TensorTraits.h b/unsupported/Eigen/CXX11/src/Tensor/TensorTraits.h
index ba09298c3..757d6cb38 100644
--- a/unsupported/Eigen/CXX11/src/Tensor/TensorTraits.h
+++ b/unsupported/Eigen/CXX11/src/Tensor/TensorTraits.h
@@ -155,6 +155,11 @@ struct eval<const TensorRef<PlainObjectType>, Eigen::Dense>
typedef const TensorRef<PlainObjectType>& type;
};
+// TODO nested<> does not exist anymore in Eigen/Core, and it thus has to be removed in favor of ref_selector.
+template<typename T, int n=1, typename PlainObject = void> struct nested
+{
+ typedef typename ref_selector<T>::type type;
+};
template <typename Scalar_, std::size_t NumIndices_, int Options_, typename IndexType_>
struct nested<Tensor<Scalar_, NumIndices_, Options_, IndexType_> >
diff --git a/unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h b/unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h
index 9e0545660..b37481cbe 100644
--- a/unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h
+++ b/unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h
@@ -400,7 +400,7 @@ template<typename Derived> struct MatrixExponentialReturnValue
Index cols() const { return m_src.cols(); }
protected:
- const typename internal::nested<Derived>::type m_src;
+ const typename internal::ref_selector<Derived>::type m_src;
};
namespace internal {
diff --git a/unsupported/Eigen/src/MatrixFunctions/MatrixFunction.h b/unsupported/Eigen/src/MatrixFunctions/MatrixFunction.h
index b68aae5e8..8f7a6f3b0 100644
--- a/unsupported/Eigen/src/MatrixFunctions/MatrixFunction.h
+++ b/unsupported/Eigen/src/MatrixFunctions/MatrixFunction.h
@@ -485,7 +485,7 @@ template<typename Derived> class MatrixFunctionReturnValue
typedef typename internal::stem_function<Scalar>::type StemFunction;
protected:
- typedef typename internal::nested<Derived>::type DerivedNested;
+ typedef typename internal::ref_selector<Derived>::type DerivedNested;
public:
diff --git a/unsupported/Eigen/src/MatrixFunctions/MatrixLogarithm.h b/unsupported/Eigen/src/MatrixFunctions/MatrixLogarithm.h
index 22bfdc4ac..463d7be0c 100644
--- a/unsupported/Eigen/src/MatrixFunctions/MatrixLogarithm.h
+++ b/unsupported/Eigen/src/MatrixFunctions/MatrixLogarithm.h
@@ -315,7 +315,7 @@ public:
typedef typename Derived::Index Index;
protected:
- typedef typename internal::nested<Derived>::type DerivedNested;
+ typedef typename internal::ref_selector<Derived>::type DerivedNested;
public:
diff --git a/unsupported/Eigen/src/MatrixFunctions/MatrixSquareRoot.h b/unsupported/Eigen/src/MatrixFunctions/MatrixSquareRoot.h
index 3a4d6eb3f..9f08c6162 100644
--- a/unsupported/Eigen/src/MatrixFunctions/MatrixSquareRoot.h
+++ b/unsupported/Eigen/src/MatrixFunctions/MatrixSquareRoot.h
@@ -320,7 +320,7 @@ template<typename Derived> class MatrixSquareRootReturnValue
{
protected:
typedef typename Derived::Index Index;
- typedef typename internal::nested<Derived>::type DerivedNested;
+ typedef typename internal::ref_selector<Derived>::type DerivedNested;
public:
/** \brief Constructor.
diff --git a/unsupported/Eigen/src/Skyline/SkylineProduct.h b/unsupported/Eigen/src/Skyline/SkylineProduct.h
index 1ddf455e2..d218a7c25 100644
--- a/unsupported/Eigen/src/Skyline/SkylineProduct.h
+++ b/unsupported/Eigen/src/Skyline/SkylineProduct.h
@@ -14,8 +14,8 @@ namespace Eigen {
template<typename Lhs, typename Rhs, int ProductMode>
struct SkylineProductReturnType {
- typedef const typename internal::nested<Lhs, Rhs::RowsAtCompileTime>::type LhsNested;
- typedef const typename internal::nested<Rhs, Lhs::RowsAtCompileTime>::type RhsNested;
+ typedef const typename internal::nested_eval<Lhs, Rhs::RowsAtCompileTime>::type LhsNested;
+ typedef const typename internal::nested_eval<Rhs, Lhs::RowsAtCompileTime>::type RhsNested;
typedef SkylineProduct<LhsNested, RhsNested, ProductMode> Type;
};
diff --git a/unsupported/Eigen/src/SparseExtra/BlockSparseMatrix.h b/unsupported/Eigen/src/SparseExtra/BlockSparseMatrix.h
index 8a7e0e57f..0e8350a7d 100644
--- a/unsupported/Eigen/src/SparseExtra/BlockSparseMatrix.h
+++ b/unsupported/Eigen/src/SparseExtra/BlockSparseMatrix.h
@@ -287,7 +287,7 @@ class BlockSparseMatrix : public SparseMatrixBase<BlockSparseMatrix<_Scalar,_Blo
typedef _Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef _StorageIndex StorageIndex;
- typedef typename internal::nested<BlockSparseMatrix<_Scalar, _BlockAtCompileTime, _Options, _StorageIndex> >::type Nested;
+ typedef typename internal::ref_selector<BlockSparseMatrix<_Scalar, _BlockAtCompileTime, _Options, _StorageIndex> >::type Nested;
enum {
Options = _Options,