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