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