aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/Core1
-rw-r--r--Eigen/src/Core/Reshape.h391
-rw-r--r--Eigen/src/Core/util/ForwardDeclarations.h1
-rw-r--r--Eigen/src/plugins/BlockMethods.h56
-rw-r--r--doc/examples/class_FixedReshape.cpp22
-rw-r--r--doc/examples/class_Reshape.cpp23
-rw-r--r--doc/snippets/MatrixBase_reshape.cpp3
-rw-r--r--doc/snippets/MatrixBase_reshape_int_int.cpp3
-rw-r--r--test/CMakeLists.txt1
-rw-r--r--test/reshape.cpp65
10 files changed, 566 insertions, 0 deletions
diff --git a/Eigen/Core b/Eigen/Core
index 9f1c63826..351c7926a 100644
--- a/Eigen/Core
+++ b/Eigen/Core
@@ -471,6 +471,7 @@ 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/Transpose.h"
diff --git a/Eigen/src/Core/Reshape.h b/Eigen/src/Core/Reshape.h
new file mode 100644
index 000000000..e4a71312c
--- /dev/null
+++ b/Eigen/src/Core/Reshape.h
@@ -0,0 +1,391 @@
+// 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) 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
+
+namespace Eigen {
+
+/** \class Reshape
+ * \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)
+ *
+ * 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
+ * most of the time this is the only way it is used.
+ *
+ * However, if you want to directly maniputate reshape expressions,
+ * for instance if you want to write a function returning such an expression, you
+ * will need to use this class.
+ *
+ * Here is an example illustrating the dynamic case:
+ * \include class_Reshape.cpp
+ * Output: \verbinclude class_Reshape.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
+ *
+ * \sa DenseBase::reshape(Index,Index), DenseBase::reshape(), class VectorReshape
+ */
+
+namespace internal {
+template<typename XprType, int ReshapeRows, int ReshapeCols>
+struct traits<Reshape<XprType, ReshapeRows, ReshapeCols> > : traits<XprType>
+{
+ typedef typename traits<XprType>::Scalar Scalar;
+ typedef typename traits<XprType>::StorageKind StorageKind;
+ typedef typename traits<XprType>::XprKind XprKind;
+ enum{
+ MatrixRows = traits<XprType>::RowsAtCompileTime,
+ MatrixCols = traits<XprType>::ColsAtCompileTime,
+ RowsAtCompileTime = ReshapeRows,
+ ColsAtCompileTime = ReshapeCols,
+ MaxRowsAtCompileTime = ReshapeRows,
+ MaxColsAtCompileTime = ReshapeCols,
+ XprTypeIsRowMajor = (int(traits<XprType>::Flags) & RowMajorBit) != 0,
+ IsRowMajor = (RowsAtCompileTime == 1 && ColsAtCompileTime != 1) ? 1
+ : (ColsAtCompileTime == 1 && RowsAtCompileTime != 1) ? 0
+ : XprTypeIsRowMajor,
+ HasSameStorageOrderAsXprType = (IsRowMajor == XprTypeIsRowMajor),
+ InnerSize = IsRowMajor ? int(ColsAtCompileTime) : int(RowsAtCompileTime),
+ InnerStrideAtCompileTime = HasSameStorageOrderAsXprType
+ ? int(inner_stride_at_compile_time<XprType>::ret)
+ : int(outer_stride_at_compile_time<XprType>::ret),
+ OuterStrideAtCompileTime = HasSameStorageOrderAsXprType
+ ? int(outer_stride_at_compile_time<XprType>::ret)
+ : int(inner_stride_at_compile_time<XprType>::ret),
+ MaskPacketAccessBit = (InnerSize == Dynamic || (InnerSize % packet_traits<Scalar>::size) == 0)
+ && (InnerStrideAtCompileTime == 1)
+ ? PacketAccessBit : 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)
+ & ~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;
+
+} // end namespace internal
+
+template<typename XprType, int ReshapeRows, int ReshapeCols, typename StorageKind> class ReshapeImpl;
+
+template<typename XprType, int ReshapeRows, int ReshapeCols> class Reshape
+ : public ReshapeImpl<XprType, ReshapeRows, ReshapeCols, typename internal::traits<XprType>::StorageKind>
+{
+ typedef ReshapeImpl<XprType, ReshapeRows, ReshapeCols, 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)
+
+ /** Fixed-size constructor
+ */
+ EIGEN_DEVICE_FUNC
+ inline Reshape(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());
+ }
+
+ /** Dynamic-size constructor
+ */
+ EIGEN_DEVICE_FUNC
+ inline Reshape(XprType& xpr,
+ Index reshapeRows, Index reshapeCols)
+ : Impl(xpr, reshapeRows, reshapeCols)
+ {
+ eigen_assert((RowsAtCompileTime==Dynamic || RowsAtCompileTime==reshapeRows)
+ && (ColsAtCompileTime==Dynamic || ColsAtCompileTime==reshapeCols));
+ eigen_assert(reshapeRows * reshapeCols == xpr.rows() * xpr.cols());
+ }
+};
+
+// The generic default implementation for dense reshape simplu forward to the internal::ReshapeImpl_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>
+{
+ typedef internal::ReshapeImpl_dense<XprType, ReshapeRows, ReshapeCols> Impl;
+ typedef typename XprType::Index Index;
+ 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)
+ : 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
+{
+ typedef Reshape<XprType, ReshapeRows, ReshapeCols> ReshapeType;
+ public:
+
+ typedef typename internal::dense_xpr_base<ReshapeType>::type Base;
+ EIGEN_DENSE_PUBLIC_INTERFACE(ReshapeType)
+ EIGEN_INHERIT_ASSIGNMENT_OPERATORS(ReshapeImpl_dense)
+
+ class InnerIterator;
+
+ /** Fixed-size constructor
+ */
+ EIGEN_DEVICE_FUNC
+ inline ReshapeImpl_dense(XprType& xpr)
+ : m_xpr(xpr), m_reshapeRows(ReshapeRows), m_reshapeCols(ReshapeCols)
+ {}
+
+ /** Dynamic-size constructor
+ */
+ EIGEN_DEVICE_FUNC
+ inline ReshapeImpl_dense(XprType& xpr,
+ Index reshapeRows, Index reshapeCols)
+ : m_xpr(xpr), m_reshapeRows(reshapeRows), m_reshapeCols(reshapeCols)
+ {}
+
+ 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
+ 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);
+ }
+
+ 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);
+ }
+
+ 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);
+
+ }
+
+ 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<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);
+ }
+
+ #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
+
+ EIGEN_DEVICE_FUNC
+ const typename internal::remove_all<typename XprType::Nested>::type& nestedExpression() const
+ {
+ return m_xpr;
+ }
+
+ protected:
+
+ 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> >
+//{
+// typedef Reshape<XprType, ReshapeRows, ReshapeCols> ReshapeType;
+// public:
+//
+// typedef MapBase<ReshapeType> Base;
+// EIGEN_DENSE_PUBLIC_INTERFACE(ReshapeType)
+// EIGEN_INHERIT_ASSIGNMENT_OPERATORS(ReshapeImpl_dense)
+//
+// /** Column or Row constructor
+// */
+// EIGEN_DEVICE_FUNC
+// inline ReshapeImpl_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()),
+// m_xpr(xpr)
+// {
+// init();
+// }
+//
+// /** Fixed-size constructor
+// */
+// EIGEN_DEVICE_FUNC
+// inline ReshapeImpl_dense(XprType& xpr)
+// : Base(internal::const_cast_ptr(&xpr.coeffRef(0, 0))), m_xpr(xpr)
+// {
+// init();
+// }
+//
+// /** Dynamic-size constructor
+// */
+// EIGEN_DEVICE_FUNC
+// inline ReshapeImpl_dense(XprType& xpr,
+// Index reshapeRows, Index reshapeCols)
+// : Base(internal::const_cast_ptr(&xpr.coeffRef(0, 0)), reshapeRows, reshapeCols),
+// m_xpr(xpr)
+// {
+// init();
+// }
+//
+// EIGEN_DEVICE_FUNC
+// const typename internal::remove_all<typename XprType::Nested>::type& nestedExpression() const
+// {
+// return m_xpr;
+// }
+//
+// EIGEN_DEVICE_FUNC
+// /** \sa MapBase::innerStride() */
+// inline Index innerStride() const
+// {
+// return internal::traits<ReshapeType>::HasSameStorageOrderAsXprType
+// ? m_xpr.innerStride()
+// : m_xpr.outerStride();
+// }
+//
+// EIGEN_DEVICE_FUNC
+// /** \sa MapBase::outerStride() */
+// inline Index outerStride() const
+// {
+// return m_outerStride;
+// }
+//
+// #ifndef __SUNPRO_CC
+// // FIXME sunstudio is not friendly with the above friend...
+// // META-FIXME there is no 'friend' keyword around here. Is this obsolete?
+// protected:
+// #endif
+//
+// #ifndef EIGEN_PARSED_BY_DOXYGEN
+// /** \internal used by allowAligned() */
+// EIGEN_DEVICE_FUNC
+// inline ReshapeImpl_dense(XprType& xpr, const Scalar* data, Index reshapeRows, Index reshapeCols)
+// : Base(data, reshapeRows, reshapeCols), m_xpr(xpr)
+// {
+// init();
+// }
+// #endif
+//
+// protected:
+// EIGEN_DEVICE_FUNC
+// void init()
+// {
+// m_outerStride = internal::traits<ReshapeType>::HasSameStorageOrderAsXprType
+// ? m_xpr.outerStride()
+// : m_xpr.innerStride();
+// }
+//
+// typename XprType::Nested m_xpr;
+// Index m_outerStride;
+//};
+
+} // end namespace internal
+
+} // end namespace Eigen
+
+#endif // EIGEN_RESHAPE_H
diff --git a/Eigen/src/Core/util/ForwardDeclarations.h b/Eigen/src/Core/util/ForwardDeclarations.h
index 1a48cff04..4e64efe24 100644
--- a/Eigen/src/Core/util/ForwardDeclarations.h
+++ b/Eigen/src/Core/util/ForwardDeclarations.h
@@ -84,6 +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 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 2d5a4e507..53401a55a 100644
--- a/Eigen/src/plugins/BlockMethods.h
+++ b/Eigen/src/plugins/BlockMethods.h
@@ -1354,3 +1354,59 @@ 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());
+}
+
diff --git a/doc/examples/class_FixedReshape.cpp b/doc/examples/class_FixedReshape.cpp
new file mode 100644
index 000000000..40156a123
--- /dev/null
+++ b/doc/examples/class_FixedReshape.cpp
@@ -0,0 +1,22 @@
+#include <Eigen/Core>
+#include <iostream>
+using namespace Eigen;
+using namespace std;
+
+template<typename Derived>
+Eigen::Reshape<Derived, 4, 2>
+reshape_helper(MatrixBase<Derived>& m)
+{
+ return Eigen::Reshape<Derived, 4, 2>(m.derived());
+}
+
+int main(int, char**)
+{
+ MatrixXd m(2, 4);
+ m << 1, 2, 3, 4,
+ 5, 6, 7, 8;
+ MatrixXd n = reshape_helper(m);
+ cout << "matrix m is:" << endl << m << endl;
+ cout << "matrix n is:" << endl << n << endl;
+ return 0;
+}
diff --git a/doc/examples/class_Reshape.cpp b/doc/examples/class_Reshape.cpp
new file mode 100644
index 000000000..b9abdc2d5
--- /dev/null
+++ b/doc/examples/class_Reshape.cpp
@@ -0,0 +1,23 @@
+#include <Eigen/Core>
+#include <iostream>
+using namespace std;
+using namespace Eigen;
+
+template<typename Derived>
+const Reshape<const Derived>
+reshape_helper(const MatrixBase<Derived>& m, int rows, int cols)
+{
+ return Reshape<const Derived>(m.derived(), rows, cols);
+}
+
+int main(int, char**)
+{
+ MatrixXd m(3, 4);
+ m << 1, 4, 7, 10,
+ 2, 5, 8, 11,
+ 3, 6, 9, 12;
+ cout << m << endl;
+ auto n = reshape_helper(m, 2, 6);
+ cout << "Matrix m is:" << endl << m << endl;
+ cout << "Matrix n is:" << endl << n << endl;
+}
diff --git a/doc/snippets/MatrixBase_reshape.cpp b/doc/snippets/MatrixBase_reshape.cpp
new file mode 100644
index 000000000..549bd1007
--- /dev/null
+++ b/doc/snippets/MatrixBase_reshape.cpp
@@ -0,0 +1,3 @@
+Matrix4i m = Matrix4i::Random();
+cout << "Here is the matrix m:" << endl << m << endl;
+cout << "Here is m.reshape<2,8>():" << endl << m.reshape<2,8>() << endl;
diff --git a/doc/snippets/MatrixBase_reshape_int_int.cpp b/doc/snippets/MatrixBase_reshape_int_int.cpp
new file mode 100644
index 000000000..1169cdb2d
--- /dev/null
+++ b/doc/snippets/MatrixBase_reshape_int_int.cpp
@@ -0,0 +1,3 @@
+Matrix4i m = Matrix4i::Random();
+cout << "Here is the matrix m:" << endl << m << endl;
+cout << "Here is m.reshape(2, 8):" << endl << m.reshape(2, 8) << endl;
diff --git a/test/CMakeLists.txt b/test/CMakeLists.txt
index 84a21b3df..315b7bece 100644
--- a/test/CMakeLists.txt
+++ b/test/CMakeLists.txt
@@ -160,6 +160,7 @@ endif()
ei_add_test(redux)
ei_add_test(visitor)
ei_add_test(block)
+ei_add_test(reshape)
ei_add_test(corners)
ei_add_test(symbolic_index)
ei_add_test(indexed_view)
diff --git a/test/reshape.cpp b/test/reshape.cpp
new file mode 100644
index 000000000..0298a2fe4
--- /dev/null
+++ b/test/reshape.cpp
@@ -0,0 +1,65 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// 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/.
+
+#include "main.h"
+
+using Eigen::Map;
+using Eigen::MatrixXi;
+
+// just test a 4x4 matrix, enumerate all combination manually,
+// so I don't have to do template-meta-programming here.
+template <typename MatType>
+void reshape_all_size(MatType m) {
+ typedef Eigen::Map<MatrixXi> MapMat;
+ // dynamic
+ VERIFY_IS_EQUAL((m.template reshape( 1, 16)), MapMat(m.data(), 1, 16));
+ VERIFY_IS_EQUAL((m.template reshape( 2, 8)), MapMat(m.data(), 2, 8));
+ VERIFY_IS_EQUAL((m.template reshape( 4, 4)), MapMat(m.data(), 4, 4));
+ VERIFY_IS_EQUAL((m.template reshape( 8, 2)), MapMat(m.data(), 8, 2));
+ VERIFY_IS_EQUAL((m.template reshape(16, 1)), MapMat(m.data(), 16, 1));
+
+ // static
+ VERIFY_IS_EQUAL((m.template reshape< 1, 16>()), MapMat(m.data(), 1, 16));
+ VERIFY_IS_EQUAL((m.template reshape< 2, 8>()), MapMat(m.data(), 2, 8));
+ VERIFY_IS_EQUAL((m.template reshape< 4, 4>()), MapMat(m.data(), 4, 4));
+ VERIFY_IS_EQUAL((m.template reshape< 8, 2>()), MapMat(m.data(), 8, 2));
+ VERIFY_IS_EQUAL((m.template reshape<16, 1>()), MapMat(m.data(), 16, 1));
+
+ // reshape chain
+ VERIFY_IS_EQUAL(
+ (m
+ .template reshape( 1, 16)
+ .template reshape< 2, 8>()
+ .template reshape(16, 1)
+ .template reshape< 8, 2>()
+ .template reshape( 2, 8)
+ .template reshape< 1, 16>()
+ .template reshape( 4, 4)
+ .template reshape<16, 1>()
+ .template reshape( 8, 2)
+ .template reshape< 4, 4>()
+ ),
+ MapMat(m.data(), 4, 4)
+ );
+}
+
+void test_reshape()
+{
+ Eigen::MatrixXi mx = Eigen::MatrixXi::Random(4, 4);
+ Eigen::Matrix4i m4 = Eigen::Matrix4i::Random(4, 4);
+
+ // test dynamic-size matrix
+ CALL_SUBTEST(reshape_all_size(mx));
+ // test static-size matrix
+ CALL_SUBTEST(reshape_all_size(m4));
+ // test dynamic-size const matrix
+ CALL_SUBTEST(reshape_all_size(static_cast<const Eigen::MatrixXi>(mx)));
+ // test static-size const matrix
+ CALL_SUBTEST(reshape_all_size(static_cast<const Eigen::Matrix4i>(m4)));
+}