Home | History | Annotate | Download | only in Core
      1 // This file is part of Eigen, a lightweight C++ template library
      2 // for linear algebra.
      3 //
      4 // Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud (at) inria.fr>
      5 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1 (at) gmail.com>
      6 //
      7 // This Source Code Form is subject to the terms of the Mozilla
      8 // Public License v. 2.0. If a copy of the MPL was not distributed
      9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
     10 
     11 #ifndef EIGEN_PARTIAL_REDUX_H
     12 #define EIGEN_PARTIAL_REDUX_H
     13 
     14 namespace Eigen {
     15 
     16 /** \class PartialReduxExpr
     17   * \ingroup Core_Module
     18   *
     19   * \brief Generic expression of a partially reduxed matrix
     20   *
     21   * \tparam MatrixType the type of the matrix we are applying the redux operation
     22   * \tparam MemberOp type of the member functor
     23   * \tparam Direction indicates the direction of the redux (#Vertical or #Horizontal)
     24   *
     25   * This class represents an expression of a partial redux operator of a matrix.
     26   * It is the return type of some VectorwiseOp functions,
     27   * and most of the time this is the only way it is used.
     28   *
     29   * \sa class VectorwiseOp
     30   */
     31 
     32 template< typename MatrixType, typename MemberOp, int Direction>
     33 class PartialReduxExpr;
     34 
     35 namespace internal {
     36 template<typename MatrixType, typename MemberOp, int Direction>
     37 struct traits<PartialReduxExpr<MatrixType, MemberOp, Direction> >
     38  : traits<MatrixType>
     39 {
     40   typedef typename MemberOp::result_type Scalar;
     41   typedef typename traits<MatrixType>::StorageKind StorageKind;
     42   typedef typename traits<MatrixType>::XprKind XprKind;
     43   typedef typename MatrixType::Scalar InputScalar;
     44   typedef typename nested<MatrixType>::type MatrixTypeNested;
     45   typedef typename remove_all<MatrixTypeNested>::type _MatrixTypeNested;
     46   enum {
     47     RowsAtCompileTime = Direction==Vertical   ? 1 : MatrixType::RowsAtCompileTime,
     48     ColsAtCompileTime = Direction==Horizontal ? 1 : MatrixType::ColsAtCompileTime,
     49     MaxRowsAtCompileTime = Direction==Vertical   ? 1 : MatrixType::MaxRowsAtCompileTime,
     50     MaxColsAtCompileTime = Direction==Horizontal ? 1 : MatrixType::MaxColsAtCompileTime,
     51     Flags0 = (unsigned int)_MatrixTypeNested::Flags & HereditaryBits,
     52     Flags = (Flags0 & ~RowMajorBit) | (RowsAtCompileTime == 1 ? RowMajorBit : 0),
     53     TraversalSize = Direction==Vertical ? RowsAtCompileTime : ColsAtCompileTime
     54   };
     55   #if EIGEN_GNUC_AT_LEAST(3,4)
     56   typedef typename MemberOp::template Cost<InputScalar,int(TraversalSize)> CostOpType;
     57   #else
     58   typedef typename MemberOp::template Cost<InputScalar,TraversalSize> CostOpType;
     59   #endif
     60   enum {
     61     CoeffReadCost = TraversalSize * traits<_MatrixTypeNested>::CoeffReadCost + int(CostOpType::value)
     62   };
     63 };
     64 }
     65 
     66 template< typename MatrixType, typename MemberOp, int Direction>
     67 class PartialReduxExpr : internal::no_assignment_operator,
     68   public internal::dense_xpr_base< PartialReduxExpr<MatrixType, MemberOp, Direction> >::type
     69 {
     70   public:
     71 
     72     typedef typename internal::dense_xpr_base<PartialReduxExpr>::type Base;
     73     EIGEN_DENSE_PUBLIC_INTERFACE(PartialReduxExpr)
     74     typedef typename internal::traits<PartialReduxExpr>::MatrixTypeNested MatrixTypeNested;
     75     typedef typename internal::traits<PartialReduxExpr>::_MatrixTypeNested _MatrixTypeNested;
     76 
     77     PartialReduxExpr(const MatrixType& mat, const MemberOp& func = MemberOp())
     78       : m_matrix(mat), m_functor(func) {}
     79 
     80     Index rows() const { return (Direction==Vertical   ? 1 : m_matrix.rows()); }
     81     Index cols() const { return (Direction==Horizontal ? 1 : m_matrix.cols()); }
     82 
     83     EIGEN_STRONG_INLINE const Scalar coeff(Index i, Index j) const
     84     {
     85       if (Direction==Vertical)
     86         return m_functor(m_matrix.col(j));
     87       else
     88         return m_functor(m_matrix.row(i));
     89     }
     90 
     91     const Scalar coeff(Index index) const
     92     {
     93       if (Direction==Vertical)
     94         return m_functor(m_matrix.col(index));
     95       else
     96         return m_functor(m_matrix.row(index));
     97     }
     98 
     99   protected:
    100     MatrixTypeNested m_matrix;
    101     const MemberOp m_functor;
    102 };
    103 
    104 #define EIGEN_MEMBER_FUNCTOR(MEMBER,COST)                               \
    105   template <typename ResultType>                                        \
    106   struct member_##MEMBER {                                           \
    107     EIGEN_EMPTY_STRUCT_CTOR(member_##MEMBER)                         \
    108     typedef ResultType result_type;                                     \
    109     template<typename Scalar, int Size> struct Cost                     \
    110     { enum { value = COST }; };                                         \
    111     template<typename XprType>                                          \
    112     EIGEN_STRONG_INLINE ResultType operator()(const XprType& mat) const \
    113     { return mat.MEMBER(); } \
    114   }
    115 
    116 namespace internal {
    117 
    118 EIGEN_MEMBER_FUNCTOR(squaredNorm, Size * NumTraits<Scalar>::MulCost + (Size-1)*NumTraits<Scalar>::AddCost);
    119 EIGEN_MEMBER_FUNCTOR(norm, (Size+5) * NumTraits<Scalar>::MulCost + (Size-1)*NumTraits<Scalar>::AddCost);
    120 EIGEN_MEMBER_FUNCTOR(stableNorm, (Size+5) * NumTraits<Scalar>::MulCost + (Size-1)*NumTraits<Scalar>::AddCost);
    121 EIGEN_MEMBER_FUNCTOR(blueNorm, (Size+5) * NumTraits<Scalar>::MulCost + (Size-1)*NumTraits<Scalar>::AddCost);
    122 EIGEN_MEMBER_FUNCTOR(hypotNorm, (Size-1) * functor_traits<scalar_hypot_op<Scalar> >::Cost );
    123 EIGEN_MEMBER_FUNCTOR(sum, (Size-1)*NumTraits<Scalar>::AddCost);
    124 EIGEN_MEMBER_FUNCTOR(mean, (Size-1)*NumTraits<Scalar>::AddCost + NumTraits<Scalar>::MulCost);
    125 EIGEN_MEMBER_FUNCTOR(minCoeff, (Size-1)*NumTraits<Scalar>::AddCost);
    126 EIGEN_MEMBER_FUNCTOR(maxCoeff, (Size-1)*NumTraits<Scalar>::AddCost);
    127 EIGEN_MEMBER_FUNCTOR(all, (Size-1)*NumTraits<Scalar>::AddCost);
    128 EIGEN_MEMBER_FUNCTOR(any, (Size-1)*NumTraits<Scalar>::AddCost);
    129 EIGEN_MEMBER_FUNCTOR(count, (Size-1)*NumTraits<Scalar>::AddCost);
    130 EIGEN_MEMBER_FUNCTOR(prod, (Size-1)*NumTraits<Scalar>::MulCost);
    131 
    132 
    133 template <typename BinaryOp, typename Scalar>
    134 struct member_redux {
    135   typedef typename result_of<
    136                      BinaryOp(Scalar)
    137                    >::type  result_type;
    138   template<typename _Scalar, int Size> struct Cost
    139   { enum { value = (Size-1) * functor_traits<BinaryOp>::Cost }; };
    140   member_redux(const BinaryOp func) : m_functor(func) {}
    141   template<typename Derived>
    142   inline result_type operator()(const DenseBase<Derived>& mat) const
    143   { return mat.redux(m_functor); }
    144   const BinaryOp m_functor;
    145 };
    146 }
    147 
    148 /** \class VectorwiseOp
    149   * \ingroup Core_Module
    150   *
    151   * \brief Pseudo expression providing partial reduction operations
    152   *
    153   * \param ExpressionType the type of the object on which to do partial reductions
    154   * \param Direction indicates the direction of the redux (#Vertical or #Horizontal)
    155   *
    156   * This class represents a pseudo expression with partial reduction features.
    157   * It is the return type of DenseBase::colwise() and DenseBase::rowwise()
    158   * and most of the time this is the only way it is used.
    159   *
    160   * Example: \include MatrixBase_colwise.cpp
    161   * Output: \verbinclude MatrixBase_colwise.out
    162   *
    163   * \sa DenseBase::colwise(), DenseBase::rowwise(), class PartialReduxExpr
    164   */
    165 template<typename ExpressionType, int Direction> class VectorwiseOp
    166 {
    167   public:
    168 
    169     typedef typename ExpressionType::Scalar Scalar;
    170     typedef typename ExpressionType::RealScalar RealScalar;
    171     typedef typename ExpressionType::Index Index;
    172     typedef typename internal::conditional<internal::must_nest_by_value<ExpressionType>::ret,
    173         ExpressionType, ExpressionType&>::type ExpressionTypeNested;
    174     typedef typename internal::remove_all<ExpressionTypeNested>::type ExpressionTypeNestedCleaned;
    175 
    176     template<template<typename _Scalar> class Functor,
    177                       typename Scalar=typename internal::traits<ExpressionType>::Scalar> struct ReturnType
    178     {
    179       typedef PartialReduxExpr<ExpressionType,
    180                                Functor<Scalar>,
    181                                Direction
    182                               > Type;
    183     };
    184 
    185     template<typename BinaryOp> struct ReduxReturnType
    186     {
    187       typedef PartialReduxExpr<ExpressionType,
    188                                internal::member_redux<BinaryOp,typename internal::traits<ExpressionType>::Scalar>,
    189                                Direction
    190                               > Type;
    191     };
    192 
    193     enum {
    194       IsVertical   = (Direction==Vertical) ? 1 : 0,
    195       IsHorizontal = (Direction==Horizontal) ? 1 : 0
    196     };
    197 
    198   protected:
    199 
    200     /** \internal
    201       * \returns the i-th subvector according to the \c Direction */
    202     typedef typename internal::conditional<Direction==Vertical,
    203                                typename ExpressionType::ColXpr,
    204                                typename ExpressionType::RowXpr>::type SubVector;
    205     SubVector subVector(Index i)
    206     {
    207       return SubVector(m_matrix.derived(),i);
    208     }
    209 
    210     /** \internal
    211       * \returns the number of subvectors in the direction \c Direction */
    212     Index subVectors() const
    213     { return Direction==Vertical?m_matrix.cols():m_matrix.rows(); }
    214 
    215     template<typename OtherDerived> struct ExtendedType {
    216       typedef Replicate<OtherDerived,
    217                         Direction==Vertical   ? 1 : ExpressionType::RowsAtCompileTime,
    218                         Direction==Horizontal ? 1 : ExpressionType::ColsAtCompileTime> Type;
    219     };
    220 
    221     /** \internal
    222       * Replicates a vector to match the size of \c *this */
    223     template<typename OtherDerived>
    224     typename ExtendedType<OtherDerived>::Type
    225     extendedTo(const DenseBase<OtherDerived>& other) const
    226     {
    227       EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(Direction==Vertical, OtherDerived::MaxColsAtCompileTime==1),
    228                           YOU_PASSED_A_ROW_VECTOR_BUT_A_COLUMN_VECTOR_WAS_EXPECTED)
    229       EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(Direction==Horizontal, OtherDerived::MaxRowsAtCompileTime==1),
    230                           YOU_PASSED_A_COLUMN_VECTOR_BUT_A_ROW_VECTOR_WAS_EXPECTED)
    231       return typename ExtendedType<OtherDerived>::Type
    232                       (other.derived(),
    233                        Direction==Vertical   ? 1 : m_matrix.rows(),
    234                        Direction==Horizontal ? 1 : m_matrix.cols());
    235     }
    236 
    237   public:
    238 
    239     inline VectorwiseOp(ExpressionType& matrix) : m_matrix(matrix) {}
    240 
    241     /** \internal */
    242     inline const ExpressionType& _expression() const { return m_matrix; }
    243 
    244     /** \returns a row or column vector expression of \c *this reduxed by \a func
    245       *
    246       * The template parameter \a BinaryOp is the type of the functor
    247       * of the custom redux operator. Note that func must be an associative operator.
    248       *
    249       * \sa class VectorwiseOp, DenseBase::colwise(), DenseBase::rowwise()
    250       */
    251     template<typename BinaryOp>
    252     const typename ReduxReturnType<BinaryOp>::Type
    253     redux(const BinaryOp& func = BinaryOp()) const
    254     { return typename ReduxReturnType<BinaryOp>::Type(_expression(), func); }
    255 
    256     /** \returns a row (or column) vector expression of the smallest coefficient
    257       * of each column (or row) of the referenced expression.
    258       *
    259       * Example: \include PartialRedux_minCoeff.cpp
    260       * Output: \verbinclude PartialRedux_minCoeff.out
    261       *
    262       * \sa DenseBase::minCoeff() */
    263     const typename ReturnType<internal::member_minCoeff>::Type minCoeff() const
    264     { return _expression(); }
    265 
    266     /** \returns a row (or column) vector expression of the largest coefficient
    267       * of each column (or row) of the referenced expression.
    268       *
    269       * Example: \include PartialRedux_maxCoeff.cpp
    270       * Output: \verbinclude PartialRedux_maxCoeff.out
    271       *
    272       * \sa DenseBase::maxCoeff() */
    273     const typename ReturnType<internal::member_maxCoeff>::Type maxCoeff() const
    274     { return _expression(); }
    275 
    276     /** \returns a row (or column) vector expression of the squared norm
    277       * of each column (or row) of the referenced expression.
    278       *
    279       * Example: \include PartialRedux_squaredNorm.cpp
    280       * Output: \verbinclude PartialRedux_squaredNorm.out
    281       *
    282       * \sa DenseBase::squaredNorm() */
    283     const typename ReturnType<internal::member_squaredNorm,RealScalar>::Type squaredNorm() const
    284     { return _expression(); }
    285 
    286     /** \returns a row (or column) vector expression of the norm
    287       * of each column (or row) of the referenced expression.
    288       *
    289       * Example: \include PartialRedux_norm.cpp
    290       * Output: \verbinclude PartialRedux_norm.out
    291       *
    292       * \sa DenseBase::norm() */
    293     const typename ReturnType<internal::member_norm,RealScalar>::Type norm() const
    294     { return _expression(); }
    295 
    296 
    297     /** \returns a row (or column) vector expression of the norm
    298       * of each column (or row) of the referenced expression, using
    299       * blue's algorithm.
    300       *
    301       * \sa DenseBase::blueNorm() */
    302     const typename ReturnType<internal::member_blueNorm,RealScalar>::Type blueNorm() const
    303     { return _expression(); }
    304 
    305 
    306     /** \returns a row (or column) vector expression of the norm
    307       * of each column (or row) of the referenced expression, avoiding
    308       * underflow and overflow.
    309       *
    310       * \sa DenseBase::stableNorm() */
    311     const typename ReturnType<internal::member_stableNorm,RealScalar>::Type stableNorm() const
    312     { return _expression(); }
    313 
    314 
    315     /** \returns a row (or column) vector expression of the norm
    316       * of each column (or row) of the referenced expression, avoiding
    317       * underflow and overflow using a concatenation of hypot() calls.
    318       *
    319       * \sa DenseBase::hypotNorm() */
    320     const typename ReturnType<internal::member_hypotNorm,RealScalar>::Type hypotNorm() const
    321     { return _expression(); }
    322 
    323     /** \returns a row (or column) vector expression of the sum
    324       * of each column (or row) of the referenced expression.
    325       *
    326       * Example: \include PartialRedux_sum.cpp
    327       * Output: \verbinclude PartialRedux_sum.out
    328       *
    329       * \sa DenseBase::sum() */
    330     const typename ReturnType<internal::member_sum>::Type sum() const
    331     { return _expression(); }
    332 
    333     /** \returns a row (or column) vector expression of the mean
    334     * of each column (or row) of the referenced expression.
    335     *
    336     * \sa DenseBase::mean() */
    337     const typename ReturnType<internal::member_mean>::Type mean() const
    338     { return _expression(); }
    339 
    340     /** \returns a row (or column) vector expression representing
    341       * whether \b all coefficients of each respective column (or row) are \c true.
    342       *
    343       * \sa DenseBase::all() */
    344     const typename ReturnType<internal::member_all>::Type all() const
    345     { return _expression(); }
    346 
    347     /** \returns a row (or column) vector expression representing
    348       * whether \b at \b least one coefficient of each respective column (or row) is \c true.
    349       *
    350       * \sa DenseBase::any() */
    351     const typename ReturnType<internal::member_any>::Type any() const
    352     { return _expression(); }
    353 
    354     /** \returns a row (or column) vector expression representing
    355       * the number of \c true coefficients of each respective column (or row).
    356       *
    357       * Example: \include PartialRedux_count.cpp
    358       * Output: \verbinclude PartialRedux_count.out
    359       *
    360       * \sa DenseBase::count() */
    361     const PartialReduxExpr<ExpressionType, internal::member_count<Index>, Direction> count() const
    362     { return _expression(); }
    363 
    364     /** \returns a row (or column) vector expression of the product
    365       * of each column (or row) of the referenced expression.
    366       *
    367       * Example: \include PartialRedux_prod.cpp
    368       * Output: \verbinclude PartialRedux_prod.out
    369       *
    370       * \sa DenseBase::prod() */
    371     const typename ReturnType<internal::member_prod>::Type prod() const
    372     { return _expression(); }
    373 
    374 
    375     /** \returns a matrix expression
    376       * where each column (or row) are reversed.
    377       *
    378       * Example: \include Vectorwise_reverse.cpp
    379       * Output: \verbinclude Vectorwise_reverse.out
    380       *
    381       * \sa DenseBase::reverse() */
    382     const Reverse<ExpressionType, Direction> reverse() const
    383     { return Reverse<ExpressionType, Direction>( _expression() ); }
    384 
    385     typedef Replicate<ExpressionType,Direction==Vertical?Dynamic:1,Direction==Horizontal?Dynamic:1> ReplicateReturnType;
    386     const ReplicateReturnType replicate(Index factor) const;
    387 
    388     /**
    389       * \return an expression of the replication of each column (or row) of \c *this
    390       *
    391       * Example: \include DirectionWise_replicate.cpp
    392       * Output: \verbinclude DirectionWise_replicate.out
    393       *
    394       * \sa VectorwiseOp::replicate(Index), DenseBase::replicate(), class Replicate
    395       */
    396     // NOTE implemented here because of sunstudio's compilation errors
    397     template<int Factor> const Replicate<ExpressionType,(IsVertical?Factor:1),(IsHorizontal?Factor:1)>
    398     replicate(Index factor = Factor) const
    399     {
    400       return Replicate<ExpressionType,Direction==Vertical?Factor:1,Direction==Horizontal?Factor:1>
    401           (_expression(),Direction==Vertical?factor:1,Direction==Horizontal?factor:1);
    402     }
    403 
    404 /////////// Artithmetic operators ///////////
    405 
    406     /** Copies the vector \a other to each subvector of \c *this */
    407     template<typename OtherDerived>
    408     ExpressionType& operator=(const DenseBase<OtherDerived>& other)
    409     {
    410       EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
    411       EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived)
    412       //eigen_assert((m_matrix.isNull()) == (other.isNull())); FIXME
    413       return const_cast<ExpressionType&>(m_matrix = extendedTo(other.derived()));
    414     }
    415 
    416     /** Adds the vector \a other to each subvector of \c *this */
    417     template<typename OtherDerived>
    418     ExpressionType& operator+=(const DenseBase<OtherDerived>& other)
    419     {
    420       EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
    421       EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived)
    422       return const_cast<ExpressionType&>(m_matrix += extendedTo(other.derived()));
    423     }
    424 
    425     /** Substracts the vector \a other to each subvector of \c *this */
    426     template<typename OtherDerived>
    427     ExpressionType& operator-=(const DenseBase<OtherDerived>& other)
    428     {
    429       EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
    430       EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived)
    431       return const_cast<ExpressionType&>(m_matrix -= extendedTo(other.derived()));
    432     }
    433 
    434     /** Multiples each subvector of \c *this by the vector \a other */
    435     template<typename OtherDerived>
    436     ExpressionType& operator*=(const DenseBase<OtherDerived>& other)
    437     {
    438       EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
    439       EIGEN_STATIC_ASSERT_ARRAYXPR(ExpressionType)
    440       EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived)
    441       m_matrix *= extendedTo(other.derived());
    442       return const_cast<ExpressionType&>(m_matrix);
    443     }
    444 
    445     /** Divides each subvector of \c *this by the vector \a other */
    446     template<typename OtherDerived>
    447     ExpressionType& operator/=(const DenseBase<OtherDerived>& other)
    448     {
    449       EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
    450       EIGEN_STATIC_ASSERT_ARRAYXPR(ExpressionType)
    451       EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived)
    452       m_matrix /= extendedTo(other.derived());
    453       return const_cast<ExpressionType&>(m_matrix);
    454     }
    455 
    456     /** Returns the expression of the sum of the vector \a other to each subvector of \c *this */
    457     template<typename OtherDerived> EIGEN_STRONG_INLINE
    458     CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
    459                   const ExpressionTypeNestedCleaned,
    460                   const typename ExtendedType<OtherDerived>::Type>
    461     operator+(const DenseBase<OtherDerived>& other) const
    462     {
    463       EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
    464       EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived)
    465       return m_matrix + extendedTo(other.derived());
    466     }
    467 
    468     /** Returns the expression of the difference between each subvector of \c *this and the vector \a other */
    469     template<typename OtherDerived>
    470     CwiseBinaryOp<internal::scalar_difference_op<Scalar>,
    471                   const ExpressionTypeNestedCleaned,
    472                   const typename ExtendedType<OtherDerived>::Type>
    473     operator-(const DenseBase<OtherDerived>& other) const
    474     {
    475       EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
    476       EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived)
    477       return m_matrix - extendedTo(other.derived());
    478     }
    479 
    480     /** Returns the expression where each subvector is the product of the vector \a other
    481       * by the corresponding subvector of \c *this */
    482     template<typename OtherDerived> EIGEN_STRONG_INLINE
    483     CwiseBinaryOp<internal::scalar_product_op<Scalar>,
    484                   const ExpressionTypeNestedCleaned,
    485                   const typename ExtendedType<OtherDerived>::Type>
    486     operator*(const DenseBase<OtherDerived>& other) const
    487     {
    488       EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
    489       EIGEN_STATIC_ASSERT_ARRAYXPR(ExpressionType)
    490       EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived)
    491       return m_matrix * extendedTo(other.derived());
    492     }
    493 
    494     /** Returns the expression where each subvector is the quotient of the corresponding
    495       * subvector of \c *this by the vector \a other */
    496     template<typename OtherDerived>
    497     CwiseBinaryOp<internal::scalar_quotient_op<Scalar>,
    498                   const ExpressionTypeNestedCleaned,
    499                   const typename ExtendedType<OtherDerived>::Type>
    500     operator/(const DenseBase<OtherDerived>& other) const
    501     {
    502       EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
    503       EIGEN_STATIC_ASSERT_ARRAYXPR(ExpressionType)
    504       EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived)
    505       return m_matrix / extendedTo(other.derived());
    506     }
    507 
    508 /////////// Geometry module ///////////
    509 
    510     #if EIGEN2_SUPPORT_STAGE > STAGE20_RESOLVE_API_CONFLICTS
    511     Homogeneous<ExpressionType,Direction> homogeneous() const;
    512     #endif
    513 
    514     typedef typename ExpressionType::PlainObject CrossReturnType;
    515     template<typename OtherDerived>
    516     const CrossReturnType cross(const MatrixBase<OtherDerived>& other) const;
    517 
    518     enum {
    519       HNormalized_Size = Direction==Vertical ? internal::traits<ExpressionType>::RowsAtCompileTime
    520                                              : internal::traits<ExpressionType>::ColsAtCompileTime,
    521       HNormalized_SizeMinusOne = HNormalized_Size==Dynamic ? Dynamic : HNormalized_Size-1
    522     };
    523     typedef Block<const ExpressionType,
    524                   Direction==Vertical   ? int(HNormalized_SizeMinusOne)
    525                                         : int(internal::traits<ExpressionType>::RowsAtCompileTime),
    526                   Direction==Horizontal ? int(HNormalized_SizeMinusOne)
    527                                         : int(internal::traits<ExpressionType>::ColsAtCompileTime)>
    528             HNormalized_Block;
    529     typedef Block<const ExpressionType,
    530                   Direction==Vertical   ? 1 : int(internal::traits<ExpressionType>::RowsAtCompileTime),
    531                   Direction==Horizontal ? 1 : int(internal::traits<ExpressionType>::ColsAtCompileTime)>
    532             HNormalized_Factors;
    533     typedef CwiseBinaryOp<internal::scalar_quotient_op<typename internal::traits<ExpressionType>::Scalar>,
    534                 const HNormalized_Block,
    535                 const Replicate<HNormalized_Factors,
    536                   Direction==Vertical   ? HNormalized_SizeMinusOne : 1,
    537                   Direction==Horizontal ? HNormalized_SizeMinusOne : 1> >
    538             HNormalizedReturnType;
    539 
    540     const HNormalizedReturnType hnormalized() const;
    541 
    542   protected:
    543     ExpressionTypeNested m_matrix;
    544 };
    545 
    546 /** \returns a VectorwiseOp wrapper of *this providing additional partial reduction operations
    547   *
    548   * Example: \include MatrixBase_colwise.cpp
    549   * Output: \verbinclude MatrixBase_colwise.out
    550   *
    551   * \sa rowwise(), class VectorwiseOp, \ref TutorialReductionsVisitorsBroadcasting
    552   */
    553 template<typename Derived>
    554 inline const typename DenseBase<Derived>::ConstColwiseReturnType
    555 DenseBase<Derived>::colwise() const
    556 {
    557   return derived();
    558 }
    559 
    560 /** \returns a writable VectorwiseOp wrapper of *this providing additional partial reduction operations
    561   *
    562   * \sa rowwise(), class VectorwiseOp, \ref TutorialReductionsVisitorsBroadcasting
    563   */
    564 template<typename Derived>
    565 inline typename DenseBase<Derived>::ColwiseReturnType
    566 DenseBase<Derived>::colwise()
    567 {
    568   return derived();
    569 }
    570 
    571 /** \returns a VectorwiseOp wrapper of *this providing additional partial reduction operations
    572   *
    573   * Example: \include MatrixBase_rowwise.cpp
    574   * Output: \verbinclude MatrixBase_rowwise.out
    575   *
    576   * \sa colwise(), class VectorwiseOp, \ref TutorialReductionsVisitorsBroadcasting
    577   */
    578 template<typename Derived>
    579 inline const typename DenseBase<Derived>::ConstRowwiseReturnType
    580 DenseBase<Derived>::rowwise() const
    581 {
    582   return derived();
    583 }
    584 
    585 /** \returns a writable VectorwiseOp wrapper of *this providing additional partial reduction operations
    586   *
    587   * \sa colwise(), class VectorwiseOp, \ref TutorialReductionsVisitorsBroadcasting
    588   */
    589 template<typename Derived>
    590 inline typename DenseBase<Derived>::RowwiseReturnType
    591 DenseBase<Derived>::rowwise()
    592 {
    593   return derived();
    594 }
    595 
    596 } // end namespace Eigen
    597 
    598 #endif // EIGEN_PARTIAL_REDUX_H
    599