Replicate.h revision c981c48f5bc9aefeffc0bcb0cc3934c2fae179dd
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_REPLICATE_H
11#define EIGEN_REPLICATE_H
12
13namespace Eigen {
14
15/**
16  * \class Replicate
17  * \ingroup Core_Module
18  *
19  * \brief Expression of the multiple replication of a matrix or vector
20  *
21  * \param MatrixType the type of the object we are replicating
22  *
23  * This class represents an expression of the multiple replication of a matrix or vector.
24  * It is the return type of DenseBase::replicate() and most of the time
25  * this is the only way it is used.
26  *
27  * \sa DenseBase::replicate()
28  */
29
30namespace internal {
31template<typename MatrixType,int RowFactor,int ColFactor>
32struct traits<Replicate<MatrixType,RowFactor,ColFactor> >
33 : traits<MatrixType>
34{
35  typedef typename MatrixType::Scalar Scalar;
36  typedef typename traits<MatrixType>::StorageKind StorageKind;
37  typedef typename traits<MatrixType>::XprKind XprKind;
38  enum {
39    Factor = (RowFactor==Dynamic || ColFactor==Dynamic) ? Dynamic : RowFactor*ColFactor
40  };
41  typedef typename nested<MatrixType,Factor>::type MatrixTypeNested;
42  typedef typename remove_reference<MatrixTypeNested>::type _MatrixTypeNested;
43  enum {
44    RowsAtCompileTime = RowFactor==Dynamic || int(MatrixType::RowsAtCompileTime)==Dynamic
45                      ? Dynamic
46                      : RowFactor * MatrixType::RowsAtCompileTime,
47    ColsAtCompileTime = ColFactor==Dynamic || int(MatrixType::ColsAtCompileTime)==Dynamic
48                      ? Dynamic
49                      : ColFactor * MatrixType::ColsAtCompileTime,
50   //FIXME we don't propagate the max sizes !!!
51    MaxRowsAtCompileTime = RowsAtCompileTime,
52    MaxColsAtCompileTime = ColsAtCompileTime,
53    IsRowMajor = MaxRowsAtCompileTime==1 && MaxColsAtCompileTime!=1 ? 1
54               : MaxColsAtCompileTime==1 && MaxRowsAtCompileTime!=1 ? 0
55               : (MatrixType::Flags & RowMajorBit) ? 1 : 0,
56    Flags = (_MatrixTypeNested::Flags & HereditaryBits & ~RowMajorBit) | (IsRowMajor ? RowMajorBit : 0),
57    CoeffReadCost = _MatrixTypeNested::CoeffReadCost
58  };
59};
60}
61
62template<typename MatrixType,int RowFactor,int ColFactor> class Replicate
63  : public internal::dense_xpr_base< Replicate<MatrixType,RowFactor,ColFactor> >::type
64{
65    typedef typename internal::traits<Replicate>::MatrixTypeNested MatrixTypeNested;
66    typedef typename internal::traits<Replicate>::_MatrixTypeNested _MatrixTypeNested;
67  public:
68
69    typedef typename internal::dense_xpr_base<Replicate>::type Base;
70    EIGEN_DENSE_PUBLIC_INTERFACE(Replicate)
71
72    template<typename OriginalMatrixType>
73    inline explicit Replicate(const OriginalMatrixType& matrix)
74      : m_matrix(matrix), m_rowFactor(RowFactor), m_colFactor(ColFactor)
75    {
76      EIGEN_STATIC_ASSERT((internal::is_same<typename internal::remove_const<MatrixType>::type,OriginalMatrixType>::value),
77                          THE_MATRIX_OR_EXPRESSION_THAT_YOU_PASSED_DOES_NOT_HAVE_THE_EXPECTED_TYPE)
78      eigen_assert(RowFactor!=Dynamic && ColFactor!=Dynamic);
79    }
80
81    template<typename OriginalMatrixType>
82    inline Replicate(const OriginalMatrixType& matrix, Index rowFactor, Index colFactor)
83      : m_matrix(matrix), m_rowFactor(rowFactor), m_colFactor(colFactor)
84    {
85      EIGEN_STATIC_ASSERT((internal::is_same<typename internal::remove_const<MatrixType>::type,OriginalMatrixType>::value),
86                          THE_MATRIX_OR_EXPRESSION_THAT_YOU_PASSED_DOES_NOT_HAVE_THE_EXPECTED_TYPE)
87    }
88
89    inline Index rows() const { return m_matrix.rows() * m_rowFactor.value(); }
90    inline Index cols() const { return m_matrix.cols() * m_colFactor.value(); }
91
92    inline Scalar coeff(Index row, Index col) const
93    {
94      // try to avoid using modulo; this is a pure optimization strategy
95      const Index actual_row  = internal::traits<MatrixType>::RowsAtCompileTime==1 ? 0
96                            : RowFactor==1 ? row
97                            : row%m_matrix.rows();
98      const Index actual_col  = internal::traits<MatrixType>::ColsAtCompileTime==1 ? 0
99                            : ColFactor==1 ? col
100                            : col%m_matrix.cols();
101
102      return m_matrix.coeff(actual_row, actual_col);
103    }
104    template<int LoadMode>
105    inline PacketScalar packet(Index row, Index col) const
106    {
107      const Index actual_row  = internal::traits<MatrixType>::RowsAtCompileTime==1 ? 0
108                            : RowFactor==1 ? row
109                            : row%m_matrix.rows();
110      const Index actual_col  = internal::traits<MatrixType>::ColsAtCompileTime==1 ? 0
111                            : ColFactor==1 ? col
112                            : col%m_matrix.cols();
113
114      return m_matrix.template packet<LoadMode>(actual_row, actual_col);
115    }
116
117    const _MatrixTypeNested& nestedExpression() const
118    {
119      return m_matrix;
120    }
121
122  protected:
123    MatrixTypeNested m_matrix;
124    const internal::variable_if_dynamic<Index, RowFactor> m_rowFactor;
125    const internal::variable_if_dynamic<Index, ColFactor> m_colFactor;
126};
127
128/**
129  * \return an expression of the replication of \c *this
130  *
131  * Example: \include MatrixBase_replicate.cpp
132  * Output: \verbinclude MatrixBase_replicate.out
133  *
134  * \sa VectorwiseOp::replicate(), DenseBase::replicate(Index,Index), class Replicate
135  */
136template<typename Derived>
137template<int RowFactor, int ColFactor>
138inline const Replicate<Derived,RowFactor,ColFactor>
139DenseBase<Derived>::replicate() const
140{
141  return Replicate<Derived,RowFactor,ColFactor>(derived());
142}
143
144/**
145  * \return an expression of the replication of \c *this
146  *
147  * Example: \include MatrixBase_replicate_int_int.cpp
148  * Output: \verbinclude MatrixBase_replicate_int_int.out
149  *
150  * \sa VectorwiseOp::replicate(), DenseBase::replicate<int,int>(), class Replicate
151  */
152template<typename Derived>
153inline const Replicate<Derived,Dynamic,Dynamic>
154DenseBase<Derived>::replicate(Index rowFactor,Index colFactor) const
155{
156  return Replicate<Derived,Dynamic,Dynamic>(derived(),rowFactor,colFactor);
157}
158
159/**
160  * \return an expression of the replication of each column (or row) of \c *this
161  *
162  * Example: \include DirectionWise_replicate_int.cpp
163  * Output: \verbinclude DirectionWise_replicate_int.out
164  *
165  * \sa VectorwiseOp::replicate(), DenseBase::replicate(), class Replicate
166  */
167template<typename ExpressionType, int Direction>
168const typename VectorwiseOp<ExpressionType,Direction>::ReplicateReturnType
169VectorwiseOp<ExpressionType,Direction>::replicate(Index factor) const
170{
171  return typename VectorwiseOp<ExpressionType,Direction>::ReplicateReturnType
172          (_expression(),Direction==Vertical?factor:1,Direction==Horizontal?factor:1);
173}
174
175} // end namespace Eigen
176
177#endif // EIGEN_REPLICATE_H
178