diff options
author | Gael Guennebaud <g.gael@free.fr> | 2017-01-29 14:57:45 +0100 |
---|---|---|
committer | Gael Guennebaud <g.gael@free.fr> | 2017-01-29 14:57:45 +0100 |
commit | 9036cda36484c4d7268b928b5976534c8ef3ce42 (patch) | |
tree | b3860ea39b51fed3ec166495e17051ece376f3c7 /Eigen | |
parent | 0e89baa5d895e40bae63c804cd3d3c568dca50f1 (diff) |
Cleanup intitial reshape implementation:
- reshape -> reshaped
- make it compatible with evaluators.
Diffstat (limited to 'Eigen')
-rw-r--r-- | Eigen/Core | 2 | ||||
-rw-r--r-- | Eigen/src/Core/DenseBase.h | 14 | ||||
-rw-r--r-- | Eigen/src/Core/Reshape.h | 441 | ||||
-rw-r--r-- | Eigen/src/Core/util/ForwardDeclarations.h | 2 | ||||
-rw-r--r-- | Eigen/src/plugins/BlockMethods.h | 56 |
5 files changed, 285 insertions, 230 deletions
diff --git a/Eigen/Core b/Eigen/Core index 351c7926a..e72e528a3 100644 --- a/Eigen/Core +++ b/Eigen/Core @@ -471,9 +471,9 @@ using std::ptrdiff_t; #include "src/Core/Map.h" #include "src/Core/Ref.h" #include "src/Core/Block.h" -#include "src/Core/Reshape.h" #include "src/Core/VectorBlock.h" #include "src/Core/IndexedView.h" +#include "src/Core/Reshape.h" #include "src/Core/Transpose.h" #include "src/Core/DiagonalMatrix.h" #include "src/Core/Diagonal.h" diff --git a/Eigen/src/Core/DenseBase.h b/Eigen/src/Core/DenseBase.h index a8229cf03..09e958727 100644 --- a/Eigen/src/Core/DenseBase.h +++ b/Eigen/src/Core/DenseBase.h @@ -557,6 +557,20 @@ template<typename Derived> class DenseBase } EIGEN_DEVICE_FUNC void reverseInPlace(); + EIGEN_DEVICE_FUNC inline + Reshaped<Derived> reshaped(Index nRows, Index nCols); + + EIGEN_DEVICE_FUNC inline + const Reshaped<const Derived> reshaped(Index nRows, Index nCols) const; + + template<int ReshapeRows, int ReshapeCols> + EIGEN_DEVICE_FUNC + inline Reshaped<Derived, ReshapeRows, ReshapeCols> reshaped(); + + template<int ReshapeRows, int ReshapeCols> + EIGEN_DEVICE_FUNC + inline const Reshaped<const Derived, ReshapeRows, ReshapeCols> reshaped() const; + #define EIGEN_CURRENT_STORAGE_BASE_CLASS Eigen::DenseBase #define EIGEN_DOC_BLOCK_ADDONS_NOT_INNER_PANEL #define EIGEN_DOC_BLOCK_ADDONS_INNER_PANEL_IF(COND) diff --git a/Eigen/src/Core/Reshape.h b/Eigen/src/Core/Reshape.h index e4a71312c..0f7d44c49 100644 --- a/Eigen/src/Core/Reshape.h +++ b/Eigen/src/Core/Reshape.h @@ -1,30 +1,30 @@ // This file is part of Eigen, a lightweight C++ template library // for linear algebra. // -// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr> -// Copyright (C) 2006-2010 Benoit Jacob <jacob.benoit.1@gmail.com> +// Copyright (C) 2008-2017 Gael Guennebaud <gael.guennebaud@inria.fr> // Copyright (C) 2014 yoco <peter.xiau@gmail.com> // // This Source Code Form is subject to the terms of the Mozilla // Public License v. 2.0. If a copy of the MPL was not distributed // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. -#ifndef EIGEN_RESHAPE_H -#define EIGEN_RESHAPE_H +#ifndef EIGEN_RESHAPED_H +#define EIGEN_RESHAPED_H namespace Eigen { -/** \class Reshape +/** \class Reshapedd * \ingroup Core_Module * * \brief Expression of a fixed-size or dynamic-size reshape * - * \param XprType the type of the expression in which we are taking a reshape - * \param ReshapeRows the number of rows of the reshape we are taking at compile time (optional) - * \param ReshapeCols the number of columns of the reshape we are taking at compile time (optional) + * \tparam XprType the type of the expression in which we are taking a reshape + * \tparam Rows the number of rows of the reshape we are taking at compile time (optional) + * \tparam Cols the number of columns of the reshape we are taking at compile time (optional) + * \tparam Order * * This class represents an expression of either a fixed-size or dynamic-size reshape. It is the return - * type of DenseBase::reshape(Index,Index) and DenseBase::reshape<int,int>() and + * type of DenseBase::reshaped(Index,Index) and DenseBase::reshape<int,int>() and * most of the time this is the only way it is used. * * However, if you want to directly maniputate reshape expressions, @@ -32,23 +32,23 @@ namespace Eigen { * will need to use this class. * * Here is an example illustrating the dynamic case: - * \include class_Reshape.cpp - * Output: \verbinclude class_Reshape.out + * \include class_Reshaped.cpp + * Output: \verbinclude class_Reshaped.out * * \note Even though this expression has dynamic size, in the case where \a XprType * has fixed size, this expression inherits a fixed maximal size which means that evaluating * it does not cause a dynamic memory allocation. * * Here is an example illustrating the fixed-size case: - * \include class_FixedReshape.cpp - * Output: \verbinclude class_FixedReshape.out + * \include class_FixedReshaped.cpp + * Output: \verbinclude class_FixedReshaped.out * - * \sa DenseBase::reshape(Index,Index), DenseBase::reshape(), class VectorReshape + * \sa DenseBase::reshaped(Index,Index), DenseBase::reshaped(), class VectorReshaped */ namespace internal { -template<typename XprType, int ReshapeRows, int ReshapeCols> -struct traits<Reshape<XprType, ReshapeRows, ReshapeCols> > : traits<XprType> +template<typename XprType, int Rows, int Cols, int Order> +struct traits<Reshaped<XprType, Rows, Cols, Order> > : traits<XprType> { typedef typename traits<XprType>::Scalar Scalar; typedef typename traits<XprType>::StorageKind StorageKind; @@ -56,10 +56,10 @@ struct traits<Reshape<XprType, ReshapeRows, ReshapeCols> > : traits<XprType> enum{ MatrixRows = traits<XprType>::RowsAtCompileTime, MatrixCols = traits<XprType>::ColsAtCompileTime, - RowsAtCompileTime = ReshapeRows, - ColsAtCompileTime = ReshapeCols, - MaxRowsAtCompileTime = ReshapeRows, - MaxColsAtCompileTime = ReshapeCols, + RowsAtCompileTime = Rows, + ColsAtCompileTime = Cols, + MaxRowsAtCompileTime = Rows, + MaxColsAtCompileTime = Cols, XprTypeIsRowMajor = (int(traits<XprType>::Flags) & RowMajorBit) != 0, IsRowMajor = (RowsAtCompileTime == 1 && ColsAtCompileTime != 1) ? 1 : (ColsAtCompileTime == 1 && RowsAtCompileTime != 1) ? 0 @@ -75,50 +75,48 @@ struct traits<Reshape<XprType, ReshapeRows, ReshapeCols> > : traits<XprType> MaskPacketAccessBit = (InnerSize == Dynamic || (InnerSize % packet_traits<Scalar>::size) == 0) && (InnerStrideAtCompileTime == 1) ? PacketAccessBit : 0, - MaskAlignedBit = ((OuterStrideAtCompileTime!=Dynamic) && (((OuterStrideAtCompileTime * int(sizeof(Scalar))) % 16) == 0)) ? AlignedBit : 0, + //MaskAlignedBit = ((OuterStrideAtCompileTime!=Dynamic) && (((OuterStrideAtCompileTime * int(sizeof(Scalar))) % 16) == 0)) ? AlignedBit : 0, FlagsLinearAccessBit = (RowsAtCompileTime == 1 || ColsAtCompileTime == 1) ? LinearAccessBit : 0, FlagsLvalueBit = is_lvalue<XprType>::value ? LvalueBit : 0, FlagsRowMajorBit = IsRowMajor ? RowMajorBit : 0, - Flags0 = traits<XprType>::Flags & ( (HereditaryBits & ~RowMajorBit) | - MaskPacketAccessBit | - MaskAlignedBit) + Flags0 = traits<XprType>::Flags & ( (HereditaryBits & ~RowMajorBit) | MaskPacketAccessBit) & ~DirectAccessBit, Flags = (Flags0 | FlagsLinearAccessBit | FlagsLvalueBit | FlagsRowMajorBit) }; }; -template<typename XprType, int ReshapeRows=Dynamic, int ReshapeCols=Dynamic, - bool HasDirectAccess = internal::has_direct_access<XprType>::ret> class ReshapeImpl_dense; +template<typename XprType, int Rows=Dynamic, int Cols=Dynamic, int Order = 0, + bool HasDirectAccess = internal::has_direct_access<XprType>::ret> class ReshapedImpl_dense; } // end namespace internal -template<typename XprType, int ReshapeRows, int ReshapeCols, typename StorageKind> class ReshapeImpl; +template<typename XprType, int Rows, int Cols, int Order, typename StorageKind> class ReshapedImpl; -template<typename XprType, int ReshapeRows, int ReshapeCols> class Reshape - : public ReshapeImpl<XprType, ReshapeRows, ReshapeCols, typename internal::traits<XprType>::StorageKind> +template<typename XprType, int Rows, int Cols, int Order> class Reshaped + : public ReshapedImpl<XprType, Rows, Cols, Order, typename internal::traits<XprType>::StorageKind> { - typedef ReshapeImpl<XprType, ReshapeRows, ReshapeCols, typename internal::traits<XprType>::StorageKind> Impl; + typedef ReshapedImpl<XprType, Rows, Cols, Order, typename internal::traits<XprType>::StorageKind> Impl; public: //typedef typename Impl::Base Base; typedef Impl Base; - EIGEN_GENERIC_PUBLIC_INTERFACE(Reshape) - EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Reshape) + EIGEN_GENERIC_PUBLIC_INTERFACE(Reshaped) + EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Reshaped) /** Fixed-size constructor */ EIGEN_DEVICE_FUNC - inline Reshape(XprType& xpr) + inline Reshaped(XprType& xpr) : Impl(xpr) { EIGEN_STATIC_ASSERT(RowsAtCompileTime!=Dynamic && ColsAtCompileTime!=Dynamic,THIS_METHOD_IS_ONLY_FOR_FIXED_SIZE) - eigen_assert(ReshapeRows * ReshapeCols == xpr.rows() * xpr.cols()); + eigen_assert(Rows * Cols == xpr.rows() * xpr.cols()); } /** Dynamic-size constructor */ EIGEN_DEVICE_FUNC - inline Reshape(XprType& xpr, + inline Reshaped(XprType& xpr, Index reshapeRows, Index reshapeCols) : Impl(xpr, reshapeRows, reshapeCols) { @@ -128,186 +126,226 @@ template<typename XprType, int ReshapeRows, int ReshapeCols> class Reshape } }; -// The generic default implementation for dense reshape simplu forward to the internal::ReshapeImpl_dense +// The generic default implementation for dense reshape simplu forward to the internal::ReshapedImpl_dense // that must be specialized for direct and non-direct access... -template<typename XprType, int ReshapeRows, int ReshapeCols> -class ReshapeImpl<XprType, ReshapeRows, ReshapeCols, Dense> - : public internal::ReshapeImpl_dense<XprType, ReshapeRows, ReshapeCols> +template<typename XprType, int Rows, int Cols, int Order> +class ReshapedImpl<XprType, Rows, Cols, Order, Dense> + : public internal::ReshapedImpl_dense<XprType, Rows, Cols, Order> { - typedef internal::ReshapeImpl_dense<XprType, ReshapeRows, ReshapeCols> Impl; - typedef typename XprType::Index Index; + typedef internal::ReshapedImpl_dense<XprType, Rows, Cols, Order> Impl; public: typedef Impl Base; - EIGEN_INHERIT_ASSIGNMENT_OPERATORS(ReshapeImpl) - EIGEN_DEVICE_FUNC inline ReshapeImpl(XprType& xpr) : Impl(xpr) {} - EIGEN_DEVICE_FUNC inline ReshapeImpl(XprType& xpr, Index reshapeRows, Index reshapeCols) + EIGEN_INHERIT_ASSIGNMENT_OPERATORS(ReshapedImpl) + EIGEN_DEVICE_FUNC inline ReshapedImpl(XprType& xpr) : Impl(xpr) {} + EIGEN_DEVICE_FUNC inline ReshapedImpl(XprType& xpr, Index reshapeRows, Index reshapeCols) : Impl(xpr, reshapeRows, reshapeCols) {} }; namespace internal { -/** \internal Internal implementation of dense Reshapes in the general case. */ -template<typename XprType, int ReshapeRows, int ReshapeCols, bool HasDirectAccess> class ReshapeImpl_dense - : public internal::dense_xpr_base<Reshape<XprType, ReshapeRows, ReshapeCols> >::type +/** \internal Internal implementation of dense Reshapeds in the general case. */ +template<typename XprType, int Rows, int Cols, int Order, bool HasDirectAccess> class ReshapedImpl_dense + : public internal::dense_xpr_base<Reshaped<XprType, Rows, Cols, Order> >::type { - typedef Reshape<XprType, ReshapeRows, ReshapeCols> ReshapeType; + typedef Reshaped<XprType, Rows, Cols, Order> ReshapedType; public: - typedef typename internal::dense_xpr_base<ReshapeType>::type Base; - EIGEN_DENSE_PUBLIC_INTERFACE(ReshapeType) - EIGEN_INHERIT_ASSIGNMENT_OPERATORS(ReshapeImpl_dense) + typedef typename internal::dense_xpr_base<ReshapedType>::type Base; + EIGEN_DENSE_PUBLIC_INTERFACE(ReshapedType) + EIGEN_INHERIT_ASSIGNMENT_OPERATORS(ReshapedImpl_dense) + + typedef typename internal::ref_selector<XprType>::non_const_type MatrixTypeNested; + typedef typename internal::remove_all<XprType>::type NestedExpression; class InnerIterator; /** Fixed-size constructor */ EIGEN_DEVICE_FUNC - inline ReshapeImpl_dense(XprType& xpr) - : m_xpr(xpr), m_reshapeRows(ReshapeRows), m_reshapeCols(ReshapeCols) + inline ReshapedImpl_dense(XprType& xpr) + : m_xpr(xpr), m_rows(Rows), m_cols(Cols) {} /** Dynamic-size constructor */ EIGEN_DEVICE_FUNC - inline ReshapeImpl_dense(XprType& xpr, - Index reshapeRows, Index reshapeCols) - : m_xpr(xpr), m_reshapeRows(reshapeRows), m_reshapeCols(reshapeCols) + inline ReshapedImpl_dense(XprType& xpr, + Index nRows, Index nCols) + : m_xpr(xpr), m_rows(nRows), m_cols(nCols) {} - EIGEN_DEVICE_FUNC inline Index rows() const { return m_reshapeRows.value(); } - EIGEN_DEVICE_FUNC inline Index cols() const { return m_reshapeCols.value(); } - - typedef std::pair<Index, Index> RowCol; - - inline RowCol index_remap(Index rowId, Index colId) const { - const Index nth_elem_idx = colId * m_reshapeRows.value() + rowId; - const Index actual_col = nth_elem_idx / m_xpr.rows(); - const Index actual_row = nth_elem_idx % m_xpr.rows(); - return RowCol(actual_row, actual_col); - } - - EIGEN_DEVICE_FUNC - inline Scalar& coeffRef(Index rowId, Index colId) - { - EIGEN_STATIC_ASSERT_LVALUE(XprType) - const RowCol row_col = index_remap(rowId, colId); - return m_xpr.const_cast_derived().coeffRef(row_col.first, row_col.second); - } - - EIGEN_DEVICE_FUNC - inline const Scalar& coeffRef(Index rowId, Index colId) const - { - const RowCol row_col = index_remap(rowId, colId); - return m_xpr.derived().coeffRef(row_col.first, row_col.second); - } - - EIGEN_DEVICE_FUNC - EIGEN_STRONG_INLINE const CoeffReturnType coeff(Index rowId, Index colId) const - { - const RowCol row_col = index_remap(rowId, colId); - return m_xpr.coeff(row_col.first, row_col.second); - } - - EIGEN_DEVICE_FUNC - inline Scalar& coeffRef(Index index) - { - EIGEN_STATIC_ASSERT_LVALUE(XprType) - const RowCol row_col = index_remap(RowsAtCompileTime == 1 ? 0 : index, - RowsAtCompileTime == 1 ? index : 0); - return m_xpr.const_cast_derived().coeffRef(row_col.first, row_col.second); + EIGEN_DEVICE_FUNC Index rows() const { return m_rows; } + EIGEN_DEVICE_FUNC Index cols() const { return m_cols; } - } + #ifdef EIGEN_PARSED_BY_DOXYGEN + /** \sa MapBase::data() */ + EIGEN_DEVICE_FUNC inline const Scalar* data() const; + EIGEN_DEVICE_FUNC inline Index innerStride() const; + EIGEN_DEVICE_FUNC inline Index outerStride() const; + #endif + /** \returns the nested expression */ EIGEN_DEVICE_FUNC - inline const Scalar& coeffRef(Index index) const - { - const RowCol row_col = index_remap(RowsAtCompileTime == 1 ? 0 : index, - RowsAtCompileTime == 1 ? index : 0); - return m_xpr.const_cast_derived().coeffRef(row_col.first, row_col.second); - } + const typename internal::remove_all<XprType>::type& + nestedExpression() const { return m_xpr; } + /** \returns the nested expression */ EIGEN_DEVICE_FUNC - inline const CoeffReturnType coeff(Index index) const - { - const RowCol row_col = index_remap(RowsAtCompileTime == 1 ? 0 : index, - RowsAtCompileTime == 1 ? index : 0); - return m_xpr.coeff(row_col.first, row_col.second); - } + typename internal::remove_reference<XprType>::type& + nestedExpression() { return m_xpr.const_cast_derived(); } - EIGEN_DEVICE_FUNC - template<int LoadMode> - inline PacketScalar packet(Index rowId, Index colId) const - { - const RowCol row_col = index_remap(rowId, colId); - return m_xpr.template packet<Unaligned>(row_col.first, row_col.second); + protected: - } + MatrixTypeNested m_xpr; + const internal::variable_if_dynamic<Index, RowsAtCompileTime> m_rows; + const internal::variable_if_dynamic<Index, ColsAtCompileTime> m_cols; +}; - template<int LoadMode> - inline void writePacket(Index rowId, Index colId, const PacketScalar& val) - { - const RowCol row_col = index_remap(rowId, colId); - m_xpr.const_cast_derived().template writePacket<Unaligned> - (row_col.first, row_col.second, val); - } - template<int LoadMode> - inline PacketScalar packet(Index index) const - { - const RowCol row_col = index_remap(RowsAtCompileTime == 1 ? 0 : index, - RowsAtCompileTime == 1 ? index : 0); - return m_xpr.template packet<Unaligned>(row_col.first, row_col.second); - } +template<typename ArgType, int Rows, int Cols, int Order> +struct unary_evaluator<Reshaped<ArgType, Rows, Cols, Order>, IndexBased> + : evaluator_base<Reshaped<ArgType, Rows, Cols, Order> > +{ + typedef Reshaped<ArgType, Rows, Cols, Order> XprType; - template<int LoadMode> - inline void writePacket(Index index, const PacketScalar& val) - { - const RowCol row_col = index_remap(RowsAtCompileTime == 1 ? 0 : index, - RowsAtCompileTime == 1 ? index : 0); - return m_xpr.template packet<Unaligned>(row_col.first, row_col.second, val); - } + enum { + CoeffReadCost = evaluator<ArgType>::CoeffReadCost /* TODO + cost of index computations */, - #ifdef EIGEN_PARSED_BY_DOXYGEN - /** \sa MapBase::data() */ - EIGEN_DEVICE_FUNC inline const Scalar* data() const; - EIGEN_DEVICE_FUNC inline Index innerStride() const; - EIGEN_DEVICE_FUNC inline Index outerStride() const; - #endif + Flags = (evaluator<ArgType>::Flags & (HereditaryBits /*| LinearAccessBit | DirectAccessBit*/)), - EIGEN_DEVICE_FUNC - const typename internal::remove_all<typename XprType::Nested>::type& nestedExpression() const - { - return m_xpr; - } + Alignment = 0 + }; - protected: + EIGEN_DEVICE_FUNC explicit unary_evaluator(const XprType& xpr) : m_argImpl(xpr.nestedExpression()), m_xpr(xpr) + { + EIGEN_INTERNAL_CHECK_COST_VALUE(CoeffReadCost); + } + + typedef typename XprType::Scalar Scalar; + typedef typename XprType::CoeffReturnType CoeffReturnType; + + typedef std::pair<Index, Index> RowCol; + + inline RowCol index_remap(Index rowId, Index colId) const { + const Index nth_elem_idx = colId * m_xpr.rows() + rowId; + const Index actual_col = nth_elem_idx / m_xpr.nestedExpression().rows(); + const Index actual_row = nth_elem_idx % m_xpr.nestedExpression().rows(); + return RowCol(actual_row, actual_col); + } + + EIGEN_DEVICE_FUNC + inline Scalar& coeffRef(Index rowId, Index colId) + { + EIGEN_STATIC_ASSERT_LVALUE(XprType) + const RowCol row_col = index_remap(rowId, colId); + return m_argImpl.coeffRef(row_col.first, row_col.second); + } + + EIGEN_DEVICE_FUNC + inline const Scalar& coeffRef(Index rowId, Index colId) const + { + const RowCol row_col = index_remap(rowId, colId); + return m_argImpl.coeffRef(row_col.first, row_col.second); + } + + EIGEN_DEVICE_FUNC + EIGEN_STRONG_INLINE const CoeffReturnType coeff(Index rowId, Index colId) const + { + const RowCol row_col = index_remap(rowId, colId); + return m_argImpl.coeff(row_col.first, row_col.second); + } + + EIGEN_DEVICE_FUNC + inline Scalar& coeffRef(Index index) + { + EIGEN_STATIC_ASSERT_LVALUE(XprType) + const RowCol row_col = index_remap(Rows == 1 ? 0 : index, + Rows == 1 ? index : 0); + return m_argImpl.coeffRef(row_col.first, row_col.second); + + } + + EIGEN_DEVICE_FUNC + inline const Scalar& coeffRef(Index index) const + { + const RowCol row_col = index_remap(Rows == 1 ? 0 : index, + Rows == 1 ? index : 0); + return m_argImpl.coeffRef(row_col.first, row_col.second); + } + + EIGEN_DEVICE_FUNC + inline const CoeffReturnType coeff(Index index) const + { + const RowCol row_col = index_remap(Rows == 1 ? 0 : index, + Rows == 1 ? index : 0); + return m_argImpl.coeff(row_col.first, row_col.second); + } +#if 0 + EIGEN_DEVICE_FUNC + template<int LoadMode> + inline PacketScalar packet(Index rowId, Index colId) const + { + const RowCol row_col = index_remap(rowId, colId); + return m_argImpl.template packet<Unaligned>(row_col.first, row_col.second); + + } + + template<int LoadMode> + EIGEN_DEVICE_FUNC + inline void writePacket(Index rowId, Index colId, const PacketScalar& val) + { + const RowCol row_col = index_remap(rowId, colId); + m_argImpl.const_cast_derived().template writePacket<Unaligned> + (row_col.first, row_col.second, val); + } + + template<int LoadMode> + EIGEN_DEVICE_FUNC + inline PacketScalar packet(Index index) const + { + const RowCol row_col = index_remap(RowsAtCompileTime == 1 ? 0 : index, + RowsAtCompileTime == 1 ? index : 0); + return m_argImpl.template packet<Unaligned>(row_col.first, row_col.second); + } + + template<int LoadMode> + EIGEN_DEVICE_FUNC + inline void writePacket(Index index, const PacketScalar& val) + { + const RowCol row_col = index_remap(RowsAtCompileTime == 1 ? 0 : index, + RowsAtCompileTime == 1 ? index : 0); + return m_argImpl.template packet<Unaligned>(row_col.first, row_col.second, val); + } +#endif +protected: + + evaluator<ArgType> m_argImpl; + const XprType& m_xpr; - const typename XprType::Nested m_xpr; - const internal::variable_if_dynamic<Index, RowsAtCompileTime> m_reshapeRows; - const internal::variable_if_dynamic<Index, ColsAtCompileTime> m_reshapeCols; }; -///** \internal Internal implementation of dense Reshapes in the direct access case.*/ -//template<typename XprType, int ReshapeRows, int ReshapeCols> -//class ReshapeImpl_dense<XprType,ReshapeRows,ReshapeCols, true> -// : public MapBase<Reshape<XprType, ReshapeRows, ReshapeCols> > + +///** \internal Internal implementation of dense Reshapeds in the direct access case.*/ +//template<typename XprType, int Rows, int Cols, int Order> +//class ReshapedImpl_dense<XprType,ReshapedRows,ReshapedCols, true> +// : public MapBase<Reshaped<XprType, Rows, Cols, Order> > //{ -// typedef Reshape<XprType, ReshapeRows, ReshapeCols> ReshapeType; +// typedef Reshaped<XprType, Rows, Cols, Order> ReshapedType; // public: // -// typedef MapBase<ReshapeType> Base; -// EIGEN_DENSE_PUBLIC_INTERFACE(ReshapeType) -// EIGEN_INHERIT_ASSIGNMENT_OPERATORS(ReshapeImpl_dense) +// typedef MapBase<ReshapedType> Base; +// EIGEN_DENSE_PUBLIC_INTERFACE(ReshapedType) +// EIGEN_INHERIT_ASSIGNMENT_OPERATORS(ReshapedImpl_dense) // // /** Column or Row constructor // */ // EIGEN_DEVICE_FUNC -// inline ReshapeImpl_dense(XprType& xpr, Index i) +// inline ReshapedImpl_dense(XprType& xpr, Index i) // : Base(internal::const_cast_ptr(&xpr.coeffRef( -// (ReshapeRows==1) && (ReshapeCols==XprType::ColsAtCompileTime) ? i : 0, -// (ReshapeRows==XprType::RowsAtCompileTime) && (ReshapeCols==1) ? i : 0)), -// ReshapeRows==1 ? 1 : xpr.rows(), -// ReshapeCols==1 ? 1 : xpr.cols()), +// (ReshapedRows==1) && (ReshapedCols==XprType::ColsAtCompileTime) ? i : 0, +// (ReshapedRows==XprType::RowsAtCompileTime) && (ReshapedCols==1) ? i : 0)), +// ReshapedRows==1 ? 1 : xpr.rows(), +// ReshapedCols==1 ? 1 : xpr.cols()), // m_xpr(xpr) // { // init(); @@ -316,7 +354,7 @@ template<typename XprType, int ReshapeRows, int ReshapeCols, bool HasDirectAcces // /** Fixed-size constructor // */ // EIGEN_DEVICE_FUNC -// inline ReshapeImpl_dense(XprType& xpr) +// inline ReshapedImpl_dense(XprType& xpr) // : Base(internal::const_cast_ptr(&xpr.coeffRef(0, 0))), m_xpr(xpr) // { // init(); @@ -325,7 +363,7 @@ template<typename XprType, int ReshapeRows, int ReshapeCols, bool HasDirectAcces // /** Dynamic-size constructor // */ // EIGEN_DEVICE_FUNC -// inline ReshapeImpl_dense(XprType& xpr, +// inline ReshapedImpl_dense(XprType& xpr, // Index reshapeRows, Index reshapeCols) // : Base(internal::const_cast_ptr(&xpr.coeffRef(0, 0)), reshapeRows, reshapeCols), // m_xpr(xpr) @@ -343,7 +381,7 @@ template<typename XprType, int ReshapeRows, int ReshapeCols, bool HasDirectAcces // /** \sa MapBase::innerStride() */ // inline Index innerStride() const // { -// return internal::traits<ReshapeType>::HasSameStorageOrderAsXprType +// return internal::traits<ReshapedType>::HasSameStorageOrderAsXprType // ? m_xpr.innerStride() // : m_xpr.outerStride(); // } @@ -364,7 +402,7 @@ template<typename XprType, int ReshapeRows, int ReshapeCols, bool HasDirectAcces // #ifndef EIGEN_PARSED_BY_DOXYGEN // /** \internal used by allowAligned() */ // EIGEN_DEVICE_FUNC -// inline ReshapeImpl_dense(XprType& xpr, const Scalar* data, Index reshapeRows, Index reshapeCols) +// inline ReshapedImpl_dense(XprType& xpr, const Scalar* data, Index reshapeRows, Index reshapeCols) // : Base(data, reshapeRows, reshapeCols), m_xpr(xpr) // { // init(); @@ -375,7 +413,7 @@ template<typename XprType, int ReshapeRows, int ReshapeCols, bool HasDirectAcces // EIGEN_DEVICE_FUNC // void init() // { -// m_outerStride = internal::traits<ReshapeType>::HasSameStorageOrderAsXprType +// m_outerStride = internal::traits<ReshapedType>::HasSameStorageOrderAsXprType // ? m_xpr.outerStride() // : m_xpr.innerStride(); // } @@ -386,6 +424,65 @@ template<typename XprType, int ReshapeRows, int ReshapeCols, bool HasDirectAcces } // end namespace internal +/** \returns a dynamic-size expression of a reshape in *this. + * + * \param reshapeRows the number of rows in the reshape + * \param reshapeCols the number of columns in the reshape + * + * Example: \include MatrixBase_reshape_int_int.cpp + * Output: \verbinclude MatrixBase_reshape_int_int.out + * + * \note Even though the returned expression has dynamic size, in the case + * when it is applied to a fixed-size matrix, it inherits a fixed maximal size, + * which means that evaluating it does not cause a dynamic memory allocation. + * + * \sa class Reshape, reshaped() + */ +template<typename Derived> +EIGEN_DEVICE_FUNC +inline Reshaped<Derived> DenseBase<Derived>::reshaped(Index reshapeRows, Index reshapeCols) +{ + return Reshaped<Derived>(derived(), reshapeRows, reshapeCols); +} + +/** This is the const version of reshaped(Index,Index). */ +template<typename Derived> +EIGEN_DEVICE_FUNC +inline const Reshaped<const Derived> DenseBase<Derived>::reshaped(Index reshapeRows, Index reshapeCols) const +{ + return Reshaped<const Derived>(derived(), reshapeRows, reshapeCols); +} + +/** \returns a fixed-size expression of a reshape in *this. + * + * The template parameters \a ReshapeRows and \a ReshapeCols are the number of + * rows and columns in the reshape. + * + * Example: \include MatrixBase_reshape.cpp + * Output: \verbinclude MatrixBase_reshape.out + * + * \note since reshape is a templated member, the keyword template has to be used + * if the matrix type is also a template parameter: \code m.template reshape<3,3>(); \endcode + * + * \sa class Reshape, reshaped(Index,Index) + */ +template<typename Derived> +template<int ReshapeRows, int ReshapeCols> +EIGEN_DEVICE_FUNC +inline Reshaped<Derived, ReshapeRows, ReshapeCols> DenseBase<Derived>::reshaped() +{ + return Reshaped<Derived, ReshapeRows, ReshapeCols>(derived()); +} + +/** This is the const version of reshape<>(Index, Index). */ +template<typename Derived> +template<int ReshapeRows, int ReshapeCols> +EIGEN_DEVICE_FUNC +inline const Reshaped<const Derived, ReshapeRows, ReshapeCols> DenseBase<Derived>::reshaped() const +{ + return Reshaped<const Derived, ReshapeRows, ReshapeCols>(derived()); +} + } // end namespace Eigen -#endif // EIGEN_RESHAPE_H +#endif // EIGEN_RESHAPED_H diff --git a/Eigen/src/Core/util/ForwardDeclarations.h b/Eigen/src/Core/util/ForwardDeclarations.h index 4e64efe24..fca8a350e 100644 --- a/Eigen/src/Core/util/ForwardDeclarations.h +++ b/Eigen/src/Core/util/ForwardDeclarations.h @@ -84,7 +84,7 @@ template<typename ExpressionType> class SwapWrapper; template<typename XprType, int BlockRows=Dynamic, int BlockCols=Dynamic, bool InnerPanel = false> class Block; template<typename XprType, typename RowIndices, typename ColIndices> class IndexedView; -template<typename XprType, int BlockRows=Dynamic, int BlockCols=Dynamic> class Reshape; +template<typename XprType, int Rows=Dynamic, int Cols=Dynamic, int Order=0> class Reshaped; template<typename MatrixType, int Size=Dynamic> class VectorBlock; template<typename MatrixType> class Transpose; diff --git a/Eigen/src/plugins/BlockMethods.h b/Eigen/src/plugins/BlockMethods.h index 53401a55a..2d5a4e507 100644 --- a/Eigen/src/plugins/BlockMethods.h +++ b/Eigen/src/plugins/BlockMethods.h @@ -1354,59 +1354,3 @@ inline typename ConstFixedSegmentReturnType<N>::Type tail(Index n = N) const EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) return typename ConstFixedSegmentReturnType<N>::Type(derived(), size() - n); } - -/** \returns a dynamic-size expression of a reshape in *this. - * - * \param reshapeRows the number of rows in the reshape - * \param reshapeCols the number of columns in the reshape - * - * Example: \include MatrixBase_reshape_int_int.cpp - * Output: \verbinclude MatrixBase_reshape_int_int.out - * - * \note Even though the returned expression has dynamic size, in the case - * when it is applied to a fixed-size matrix, it inherits a fixed maximal size, - * which means that evaluating it does not cause a dynamic memory allocation. - * - * \sa class Reshape, reshape() - */ -EIGEN_DEVICE_FUNC -inline Reshape<Derived> reshape(Index reshapeRows, Index reshapeCols) -{ - return Reshape<Derived>(derived(), reshapeRows, reshapeCols); -} - -/** This is the const version of reshape(Index,Index). */ -EIGEN_DEVICE_FUNC -inline const Reshape<const Derived> reshape(Index reshapeRows, Index reshapeCols) const -{ - return Reshape<const Derived>(derived(), reshapeRows, reshapeCols); -} - -/** \returns a fixed-size expression of a reshape in *this. - * - * The template parameters \a ReshapeRows and \a ReshapeCols are the number of - * rows and columns in the reshape. - * - * Example: \include MatrixBase_reshape.cpp - * Output: \verbinclude MatrixBase_reshape.out - * - * \note since reshape is a templated member, the keyword template has to be used - * if the matrix type is also a template parameter: \code m.template reshape<3,3>(); \endcode - * - * \sa class Reshape, reshape(Index,Index) - */ -template<int ReshapeRows, int ReshapeCols> -EIGEN_DEVICE_FUNC -inline Reshape<Derived, ReshapeRows, ReshapeCols> reshape() -{ - return Reshape<Derived, ReshapeRows, ReshapeCols>(derived()); -} - -/** This is the const version of reshape<>(Index, Index). */ -template<int ReshapeRows, int ReshapeCols> -EIGEN_DEVICE_FUNC -inline const Reshape<const Derived, ReshapeRows, ReshapeCols> reshape() const -{ - return Reshape<const Derived, ReshapeRows, ReshapeCols>(derived()); -} - |