Important changes to repositories hosted on mbed.com
Mbed hosted mercurial repositories are deprecated and are due to be permanently deleted in July 2026.
To keep a copy of this software download the repository Zip archive or clone locally using Mercurial.
It is also possible to export all your personal repositories from the account settings page.
VectorwiseOp.h
00001 // This file is part of Eigen, a lightweight C++ template library 00002 // for linear algebra. 00003 // 00004 // Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud@inria.fr> 00005 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com> 00006 // 00007 // This Source Code Form is subject to the terms of the Mozilla 00008 // Public License v. 2.0. If a copy of the MPL was not distributed 00009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 00010 00011 #ifndef EIGEN_PARTIAL_REDUX_H 00012 #define EIGEN_PARTIAL_REDUX_H 00013 00014 namespace Eigen { 00015 00016 /** \class PartialReduxExpr 00017 * \ingroup Core_Module 00018 * 00019 * \brief Generic expression of a partially reduxed matrix 00020 * 00021 * \tparam MatrixType the type of the matrix we are applying the redux operation 00022 * \tparam MemberOp type of the member functor 00023 * \tparam Direction indicates the direction of the redux (#Vertical or #Horizontal) 00024 * 00025 * This class represents an expression of a partial redux operator of a matrix. 00026 * It is the return type of some VectorwiseOp functions, 00027 * and most of the time this is the only way it is used. 00028 * 00029 * \sa class VectorwiseOp 00030 */ 00031 00032 template< typename MatrixType, typename MemberOp, int Direction> 00033 class PartialReduxExpr; 00034 00035 namespace internal { 00036 template<typename MatrixType, typename MemberOp, int Direction> 00037 struct traits<PartialReduxExpr<MatrixType, MemberOp, Direction> > 00038 : traits<MatrixType> 00039 { 00040 typedef typename MemberOp::result_type Scalar; 00041 typedef typename traits<MatrixType>::StorageKind StorageKind; 00042 typedef typename traits<MatrixType>::XprKind XprKind; 00043 typedef typename MatrixType::Scalar InputScalar; 00044 typedef typename nested<MatrixType>::type MatrixTypeNested; 00045 typedef typename remove_all<MatrixTypeNested>::type _MatrixTypeNested; 00046 enum { 00047 RowsAtCompileTime = Direction==Vertical ? 1 : MatrixType::RowsAtCompileTime, 00048 ColsAtCompileTime = Direction==Horizontal ? 1 : MatrixType::ColsAtCompileTime, 00049 MaxRowsAtCompileTime = Direction==Vertical ? 1 : MatrixType::MaxRowsAtCompileTime, 00050 MaxColsAtCompileTime = Direction==Horizontal ? 1 : MatrixType::MaxColsAtCompileTime, 00051 Flags0 = (unsigned int)_MatrixTypeNested::Flags & HereditaryBits, 00052 Flags = (Flags0 & ~RowMajorBit) | (RowsAtCompileTime == 1 ? RowMajorBit : 0), 00053 TraversalSize = Direction==Vertical ? MatrixType::RowsAtCompileTime : MatrixType::ColsAtCompileTime 00054 }; 00055 #if EIGEN_GNUC_AT_LEAST(3,4) 00056 typedef typename MemberOp::template Cost<InputScalar,int(TraversalSize)> CostOpType; 00057 #else 00058 typedef typename MemberOp::template Cost<InputScalar,TraversalSize> CostOpType; 00059 #endif 00060 enum { 00061 CoeffReadCost = TraversalSize==Dynamic ? Dynamic 00062 : TraversalSize * traits<_MatrixTypeNested>::CoeffReadCost + int(CostOpType::value) 00063 }; 00064 }; 00065 } 00066 00067 template< typename MatrixType, typename MemberOp, int Direction> 00068 class PartialReduxExpr : internal::no_assignment_operator, 00069 public internal::dense_xpr_base< PartialReduxExpr<MatrixType, MemberOp, Direction> >::type 00070 { 00071 public: 00072 00073 typedef typename internal::dense_xpr_base<PartialReduxExpr>::type Base; 00074 EIGEN_DENSE_PUBLIC_INTERFACE(PartialReduxExpr) 00075 typedef typename internal::traits<PartialReduxExpr>::MatrixTypeNested MatrixTypeNested; 00076 typedef typename internal::traits<PartialReduxExpr>::_MatrixTypeNested _MatrixTypeNested; 00077 00078 PartialReduxExpr(const MatrixType& mat, const MemberOp& func = MemberOp()) 00079 : m_matrix(mat), m_functor(func) {} 00080 00081 Index rows() const { return (Direction==Vertical ? 1 : m_matrix.rows()); } 00082 Index cols() const { return (Direction==Horizontal ? 1 : m_matrix.cols()); } 00083 00084 EIGEN_STRONG_INLINE const Scalar coeff(Index i, Index j) const 00085 { 00086 if (Direction==Vertical) 00087 return m_functor(m_matrix.col(j)); 00088 else 00089 return m_functor(m_matrix.row(i)); 00090 } 00091 00092 const Scalar coeff(Index index) const 00093 { 00094 if (Direction==Vertical) 00095 return m_functor(m_matrix.col(index)); 00096 else 00097 return m_functor(m_matrix.row(index)); 00098 } 00099 00100 protected: 00101 MatrixTypeNested m_matrix; 00102 const MemberOp m_functor; 00103 }; 00104 00105 #define EIGEN_MEMBER_FUNCTOR(MEMBER,COST) \ 00106 template <typename ResultType> \ 00107 struct member_##MEMBER { \ 00108 EIGEN_EMPTY_STRUCT_CTOR(member_##MEMBER) \ 00109 typedef ResultType result_type; \ 00110 template<typename Scalar, int Size> struct Cost \ 00111 { enum { value = COST }; }; \ 00112 template<typename XprType> \ 00113 EIGEN_STRONG_INLINE ResultType operator()(const XprType& mat) const \ 00114 { return mat.MEMBER(); } \ 00115 } 00116 00117 namespace internal { 00118 00119 EIGEN_MEMBER_FUNCTOR(squaredNorm, Size * NumTraits<Scalar>::MulCost + (Size-1)*NumTraits<Scalar>::AddCost); 00120 EIGEN_MEMBER_FUNCTOR(norm, (Size+5) * NumTraits<Scalar>::MulCost + (Size-1)*NumTraits<Scalar>::AddCost); 00121 EIGEN_MEMBER_FUNCTOR(stableNorm, (Size+5) * NumTraits<Scalar>::MulCost + (Size-1)*NumTraits<Scalar>::AddCost); 00122 EIGEN_MEMBER_FUNCTOR(blueNorm, (Size+5) * NumTraits<Scalar>::MulCost + (Size-1)*NumTraits<Scalar>::AddCost); 00123 EIGEN_MEMBER_FUNCTOR(hypotNorm, (Size-1) * functor_traits<scalar_hypot_op<Scalar> >::Cost ); 00124 EIGEN_MEMBER_FUNCTOR(sum, (Size-1)*NumTraits<Scalar>::AddCost); 00125 EIGEN_MEMBER_FUNCTOR(mean, (Size-1)*NumTraits<Scalar>::AddCost + NumTraits<Scalar>::MulCost); 00126 EIGEN_MEMBER_FUNCTOR(minCoeff, (Size-1)*NumTraits<Scalar>::AddCost); 00127 EIGEN_MEMBER_FUNCTOR(maxCoeff, (Size-1)*NumTraits<Scalar>::AddCost); 00128 EIGEN_MEMBER_FUNCTOR(all, (Size-1)*NumTraits<Scalar>::AddCost); 00129 EIGEN_MEMBER_FUNCTOR(any, (Size-1)*NumTraits<Scalar>::AddCost); 00130 EIGEN_MEMBER_FUNCTOR(count, (Size-1)*NumTraits<Scalar>::AddCost); 00131 EIGEN_MEMBER_FUNCTOR(prod, (Size-1)*NumTraits<Scalar>::MulCost); 00132 00133 00134 template <typename BinaryOp, typename Scalar> 00135 struct member_redux { 00136 typedef typename result_of< 00137 BinaryOp(Scalar) 00138 >::type result_type; 00139 template<typename _Scalar, int Size> struct Cost 00140 { enum { value = (Size-1) * functor_traits<BinaryOp>::Cost }; }; 00141 member_redux(const BinaryOp func) : m_functor(func) {} 00142 template<typename Derived> 00143 inline result_type operator()(const DenseBase<Derived> & mat) const 00144 { return mat.redux(m_functor); } 00145 const BinaryOp m_functor; 00146 }; 00147 } 00148 00149 /** \class VectorwiseOp 00150 * \ingroup Core_Module 00151 * 00152 * \brief Pseudo expression providing partial reduction operations 00153 * 00154 * \param ExpressionType the type of the object on which to do partial reductions 00155 * \param Direction indicates the direction of the redux (#Vertical or #Horizontal) 00156 * 00157 * This class represents a pseudo expression with partial reduction features. 00158 * It is the return type of DenseBase::colwise() and DenseBase::rowwise() 00159 * and most of the time this is the only way it is used. 00160 * 00161 * Example: \include MatrixBase_colwise.cpp 00162 * Output: \verbinclude MatrixBase_colwise.out 00163 * 00164 * \sa DenseBase::colwise(), DenseBase::rowwise(), class PartialReduxExpr 00165 */ 00166 template<typename ExpressionType, int Direction> class VectorwiseOp 00167 { 00168 public: 00169 00170 typedef typename ExpressionType::Scalar Scalar; 00171 typedef typename ExpressionType::RealScalar RealScalar; 00172 typedef typename ExpressionType::Index Index; 00173 typedef typename internal::conditional<internal::must_nest_by_value<ExpressionType>::ret, 00174 ExpressionType, ExpressionType&>::type ExpressionTypeNested; 00175 typedef typename internal::remove_all<ExpressionTypeNested>::type ExpressionTypeNestedCleaned; 00176 00177 template<template<typename _Scalar> class Functor, 00178 typename Scalar=typename internal::traits<ExpressionType>::Scalar> struct ReturnType 00179 { 00180 typedef PartialReduxExpr<ExpressionType, 00181 Functor<Scalar>, 00182 Direction 00183 > Type; 00184 }; 00185 00186 template<typename BinaryOp> struct ReduxReturnType 00187 { 00188 typedef PartialReduxExpr<ExpressionType, 00189 internal::member_redux<BinaryOp,typename internal::traits<ExpressionType>::Scalar>, 00190 Direction 00191 > Type; 00192 }; 00193 00194 enum { 00195 IsVertical = (Direction==Vertical) ? 1 : 0, 00196 IsHorizontal = (Direction==Horizontal) ? 1 : 0 00197 }; 00198 00199 protected: 00200 00201 /** \internal 00202 * \returns the i-th subvector according to the \c Direction */ 00203 typedef typename internal::conditional<Direction==Vertical, 00204 typename ExpressionType::ColXpr, 00205 typename ExpressionType::RowXpr>::type SubVector; 00206 SubVector subVector(Index i) 00207 { 00208 return SubVector(m_matrix.derived(),i); 00209 } 00210 00211 /** \internal 00212 * \returns the number of subvectors in the direction \c Direction */ 00213 Index subVectors() const 00214 { return Direction==Vertical?m_matrix.cols():m_matrix.rows(); } 00215 00216 template<typename OtherDerived> struct ExtendedType { 00217 typedef Replicate<OtherDerived, 00218 Direction==Vertical ? 1 : ExpressionType::RowsAtCompileTime, 00219 Direction==Horizontal ? 1 : ExpressionType::ColsAtCompileTime> Type; 00220 }; 00221 00222 /** \internal 00223 * Replicates a vector to match the size of \c *this */ 00224 template<typename OtherDerived> 00225 typename ExtendedType<OtherDerived>::Type 00226 extendedTo(const DenseBase<OtherDerived>& other) const 00227 { 00228 EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(Direction==Vertical, OtherDerived::MaxColsAtCompileTime==1), 00229 YOU_PASSED_A_ROW_VECTOR_BUT_A_COLUMN_VECTOR_WAS_EXPECTED) 00230 EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(Direction==Horizontal, OtherDerived::MaxRowsAtCompileTime==1), 00231 YOU_PASSED_A_COLUMN_VECTOR_BUT_A_ROW_VECTOR_WAS_EXPECTED) 00232 return typename ExtendedType<OtherDerived>::Type 00233 (other.derived(), 00234 Direction==Vertical ? 1 : m_matrix.rows(), 00235 Direction==Horizontal ? 1 : m_matrix.cols()); 00236 } 00237 00238 template<typename OtherDerived> struct OppositeExtendedType { 00239 typedef Replicate<OtherDerived, 00240 Direction==Horizontal ? 1 : ExpressionType::RowsAtCompileTime, 00241 Direction==Vertical ? 1 : ExpressionType::ColsAtCompileTime> Type; 00242 }; 00243 00244 /** \internal 00245 * Replicates a vector in the opposite direction to match the size of \c *this */ 00246 template<typename OtherDerived> 00247 typename OppositeExtendedType<OtherDerived>::Type 00248 extendedToOpposite(const DenseBase<OtherDerived>& other) const 00249 { 00250 EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(Direction==Horizontal, OtherDerived::MaxColsAtCompileTime==1), 00251 YOU_PASSED_A_ROW_VECTOR_BUT_A_COLUMN_VECTOR_WAS_EXPECTED) 00252 EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(Direction==Vertical, OtherDerived::MaxRowsAtCompileTime==1), 00253 YOU_PASSED_A_COLUMN_VECTOR_BUT_A_ROW_VECTOR_WAS_EXPECTED) 00254 return typename OppositeExtendedType<OtherDerived>::Type 00255 (other.derived(), 00256 Direction==Horizontal ? 1 : m_matrix.rows(), 00257 Direction==Vertical ? 1 : m_matrix.cols()); 00258 } 00259 00260 public: 00261 00262 inline VectorwiseOp(ExpressionType& matrix) : m_matrix(matrix) {} 00263 00264 /** \internal */ 00265 inline const ExpressionType& _expression() const { return m_matrix; } 00266 00267 /** \returns a row or column vector expression of \c *this reduxed by \a func 00268 * 00269 * The template parameter \a BinaryOp is the type of the functor 00270 * of the custom redux operator. Note that func must be an associative operator. 00271 * 00272 * \sa class VectorwiseOp, DenseBase::colwise(), DenseBase::rowwise() 00273 */ 00274 template<typename BinaryOp> 00275 const typename ReduxReturnType<BinaryOp>::Type 00276 redux(const BinaryOp& func = BinaryOp()) const 00277 { return typename ReduxReturnType<BinaryOp>::Type(_expression(), func); } 00278 00279 /** \returns a row (or column) vector expression of the smallest coefficient 00280 * of each column (or row) of the referenced expression. 00281 * 00282 * \warning the result is undefined if \c *this contains NaN. 00283 * 00284 * Example: \include PartialRedux_minCoeff.cpp 00285 * Output: \verbinclude PartialRedux_minCoeff.out 00286 * 00287 * \sa DenseBase::minCoeff() */ 00288 const typename ReturnType<internal::member_minCoeff>::Type minCoeff() const 00289 { return _expression(); } 00290 00291 /** \returns a row (or column) vector expression of the largest coefficient 00292 * of each column (or row) of the referenced expression. 00293 * 00294 * \warning the result is undefined if \c *this contains NaN. 00295 * 00296 * Example: \include PartialRedux_maxCoeff.cpp 00297 * Output: \verbinclude PartialRedux_maxCoeff.out 00298 * 00299 * \sa DenseBase::maxCoeff() */ 00300 const typename ReturnType<internal::member_maxCoeff>::Type maxCoeff() const 00301 { return _expression(); } 00302 00303 /** \returns a row (or column) vector expression of the squared norm 00304 * of each column (or row) of the referenced expression. 00305 * 00306 * Example: \include PartialRedux_squaredNorm.cpp 00307 * Output: \verbinclude PartialRedux_squaredNorm.out 00308 * 00309 * \sa DenseBase::squaredNorm() */ 00310 const typename ReturnType<internal::member_squaredNorm,RealScalar>::Type squaredNorm() const 00311 { return _expression(); } 00312 00313 /** \returns a row (or column) vector expression of the norm 00314 * of each column (or row) of the referenced expression. 00315 * 00316 * Example: \include PartialRedux_norm.cpp 00317 * Output: \verbinclude PartialRedux_norm.out 00318 * 00319 * \sa DenseBase::norm() */ 00320 const typename ReturnType<internal::member_norm,RealScalar>::Type norm() const 00321 { return _expression(); } 00322 00323 00324 /** \returns a row (or column) vector expression of the norm 00325 * of each column (or row) of the referenced expression, using 00326 * blue's algorithm. 00327 * 00328 * \sa DenseBase::blueNorm() */ 00329 const typename ReturnType<internal::member_blueNorm,RealScalar>::Type blueNorm() const 00330 { return _expression(); } 00331 00332 00333 /** \returns a row (or column) vector expression of the norm 00334 * of each column (or row) of the referenced expression, avoiding 00335 * underflow and overflow. 00336 * 00337 * \sa DenseBase::stableNorm() */ 00338 const typename ReturnType<internal::member_stableNorm,RealScalar>::Type stableNorm() const 00339 { return _expression(); } 00340 00341 00342 /** \returns a row (or column) vector expression of the norm 00343 * of each column (or row) of the referenced expression, avoiding 00344 * underflow and overflow using a concatenation of hypot() calls. 00345 * 00346 * \sa DenseBase::hypotNorm() */ 00347 const typename ReturnType<internal::member_hypotNorm,RealScalar>::Type hypotNorm() const 00348 { return _expression(); } 00349 00350 /** \returns a row (or column) vector expression of the sum 00351 * of each column (or row) of the referenced expression. 00352 * 00353 * Example: \include PartialRedux_sum.cpp 00354 * Output: \verbinclude PartialRedux_sum.out 00355 * 00356 * \sa DenseBase::sum() */ 00357 const typename ReturnType<internal::member_sum>::Type sum() const 00358 { return _expression(); } 00359 00360 /** \returns a row (or column) vector expression of the mean 00361 * of each column (or row) of the referenced expression. 00362 * 00363 * \sa DenseBase::mean() */ 00364 const typename ReturnType<internal::member_mean>::Type mean() const 00365 { return _expression(); } 00366 00367 /** \returns a row (or column) vector expression representing 00368 * whether \b all coefficients of each respective column (or row) are \c true. 00369 * 00370 * \sa DenseBase::all() */ 00371 const typename ReturnType<internal::member_all>::Type all() const 00372 { return _expression(); } 00373 00374 /** \returns a row (or column) vector expression representing 00375 * whether \b at \b least one coefficient of each respective column (or row) is \c true. 00376 * 00377 * \sa DenseBase::any() */ 00378 const typename ReturnType<internal::member_any>::Type any() const 00379 { return _expression(); } 00380 00381 /** \returns a row (or column) vector expression representing 00382 * the number of \c true coefficients of each respective column (or row). 00383 * 00384 * Example: \include PartialRedux_count.cpp 00385 * Output: \verbinclude PartialRedux_count.out 00386 * 00387 * \sa DenseBase::count() */ 00388 const PartialReduxExpr<ExpressionType, internal::member_count<Index>, Direction> count() const 00389 { return _expression(); } 00390 00391 /** \returns a row (or column) vector expression of the product 00392 * of each column (or row) of the referenced expression. 00393 * 00394 * Example: \include PartialRedux_prod.cpp 00395 * Output: \verbinclude PartialRedux_prod.out 00396 * 00397 * \sa DenseBase::prod() */ 00398 const typename ReturnType<internal::member_prod>::Type prod() const 00399 { return _expression(); } 00400 00401 00402 /** \returns a matrix expression 00403 * where each column (or row) are reversed. 00404 * 00405 * Example: \include Vectorwise_reverse.cpp 00406 * Output: \verbinclude Vectorwise_reverse.out 00407 * 00408 * \sa DenseBase::reverse() */ 00409 const Reverse<ExpressionType, Direction> reverse() const 00410 { return Reverse<ExpressionType, Direction>( _expression() ); } 00411 00412 typedef Replicate<ExpressionType,Direction==Vertical?Dynamic:1,Direction==Horizontal?Dynamic:1> ReplicateReturnType; 00413 const ReplicateReturnType replicate(Index factor) const; 00414 00415 /** 00416 * \return an expression of the replication of each column (or row) of \c *this 00417 * 00418 * Example: \include DirectionWise_replicate.cpp 00419 * Output: \verbinclude DirectionWise_replicate.out 00420 * 00421 * \sa VectorwiseOp::replicate(Index), DenseBase::replicate(), class Replicate 00422 */ 00423 // NOTE implemented here because of sunstudio's compilation errors 00424 template<int Factor> const Replicate<ExpressionType,(IsVertical?Factor:1),(IsHorizontal?Factor:1)> 00425 replicate(Index factor = Factor) const 00426 { 00427 return Replicate<ExpressionType,Direction==Vertical?Factor:1,Direction==Horizontal?Factor:1> 00428 (_expression(),Direction==Vertical?factor:1,Direction==Horizontal?factor:1); 00429 } 00430 00431 /////////// Artithmetic operators /////////// 00432 00433 /** Copies the vector \a other to each subvector of \c *this */ 00434 template<typename OtherDerived> 00435 ExpressionType& operator=(const DenseBase<OtherDerived>& other) 00436 { 00437 EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived) 00438 EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived) 00439 //eigen_assert((m_matrix.isNull()) == (other.isNull())); FIXME 00440 return const_cast<ExpressionType&>(m_matrix = extendedTo(other.derived())); 00441 } 00442 00443 /** Adds the vector \a other to each subvector of \c *this */ 00444 template<typename OtherDerived> 00445 ExpressionType& operator+=(const DenseBase<OtherDerived>& other) 00446 { 00447 EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived) 00448 EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived) 00449 return const_cast<ExpressionType&>(m_matrix += extendedTo(other.derived())); 00450 } 00451 00452 /** Substracts the vector \a other to each subvector of \c *this */ 00453 template<typename OtherDerived> 00454 ExpressionType& operator-=(const DenseBase<OtherDerived>& other) 00455 { 00456 EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived) 00457 EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived) 00458 return const_cast<ExpressionType&>(m_matrix -= extendedTo(other.derived())); 00459 } 00460 00461 /** Multiples each subvector of \c *this by the vector \a other */ 00462 template<typename OtherDerived> 00463 ExpressionType& operator*=(const DenseBase<OtherDerived>& other) 00464 { 00465 EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived) 00466 EIGEN_STATIC_ASSERT_ARRAYXPR(ExpressionType) 00467 EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived) 00468 m_matrix *= extendedTo(other.derived()); 00469 return const_cast<ExpressionType&>(m_matrix); 00470 } 00471 00472 /** Divides each subvector of \c *this by the vector \a other */ 00473 template<typename OtherDerived> 00474 ExpressionType& operator/=(const DenseBase<OtherDerived>& other) 00475 { 00476 EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived) 00477 EIGEN_STATIC_ASSERT_ARRAYXPR(ExpressionType) 00478 EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived) 00479 m_matrix /= extendedTo(other.derived()); 00480 return const_cast<ExpressionType&>(m_matrix); 00481 } 00482 00483 /** Returns the expression of the sum of the vector \a other to each subvector of \c *this */ 00484 template<typename OtherDerived> EIGEN_STRONG_INLINE 00485 CwiseBinaryOp<internal::scalar_sum_op<Scalar>, 00486 const ExpressionTypeNestedCleaned, 00487 const typename ExtendedType<OtherDerived>::Type> 00488 operator+(const DenseBase<OtherDerived>& other) const 00489 { 00490 EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived) 00491 EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived) 00492 return m_matrix + extendedTo(other.derived()); 00493 } 00494 00495 /** Returns the expression of the difference between each subvector of \c *this and the vector \a other */ 00496 template<typename OtherDerived> 00497 CwiseBinaryOp<internal::scalar_difference_op<Scalar>, 00498 const ExpressionTypeNestedCleaned, 00499 const typename ExtendedType<OtherDerived>::Type> 00500 operator-(const DenseBase<OtherDerived>& other) const 00501 { 00502 EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived) 00503 EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived) 00504 return m_matrix - extendedTo(other.derived()); 00505 } 00506 00507 /** Returns the expression where each subvector is the product of the vector \a other 00508 * by the corresponding subvector of \c *this */ 00509 template<typename OtherDerived> EIGEN_STRONG_INLINE 00510 CwiseBinaryOp<internal::scalar_product_op<Scalar>, 00511 const ExpressionTypeNestedCleaned, 00512 const typename ExtendedType<OtherDerived>::Type> 00513 operator* (const DenseBase<OtherDerived>& other) const 00514 { 00515 EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived) 00516 EIGEN_STATIC_ASSERT_ARRAYXPR(ExpressionType) 00517 EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived) 00518 return m_matrix * extendedTo(other.derived()); 00519 } 00520 00521 /** Returns the expression where each subvector is the quotient of the corresponding 00522 * subvector of \c *this by the vector \a other */ 00523 template<typename OtherDerived> 00524 CwiseBinaryOp<internal::scalar_quotient_op<Scalar>, 00525 const ExpressionTypeNestedCleaned, 00526 const typename ExtendedType<OtherDerived>::Type> 00527 operator/(const DenseBase<OtherDerived>& other) const 00528 { 00529 EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived) 00530 EIGEN_STATIC_ASSERT_ARRAYXPR(ExpressionType) 00531 EIGEN_STATIC_ASSERT_SAME_XPR_KIND(ExpressionType, OtherDerived) 00532 return m_matrix / extendedTo(other.derived()); 00533 } 00534 00535 /** \returns an expression where each column of row of the referenced matrix are normalized. 00536 * The referenced matrix is \b not modified. 00537 * \sa MatrixBase::normalized(), normalize() 00538 */ 00539 CwiseBinaryOp<internal::scalar_quotient_op<Scalar>, 00540 const ExpressionTypeNestedCleaned, 00541 const typename OppositeExtendedType<typename ReturnType<internal::member_norm,RealScalar>::Type>::Type> 00542 normalized () const { return m_matrix.cwiseQuotient(extendedToOpposite(this->norm())); } 00543 00544 00545 /** Normalize in-place each row or columns of the referenced matrix. 00546 * \sa MatrixBase::normalize(), normalized() 00547 */ 00548 void normalize() { 00549 m_matrix = this->normalized(); 00550 } 00551 00552 /////////// Geometry module /////////// 00553 00554 #if EIGEN2_SUPPORT_STAGE > STAGE20_RESOLVE_API_CONFLICTS 00555 Homogeneous<ExpressionType,Direction> homogeneous() const; 00556 #endif 00557 00558 typedef typename ExpressionType::PlainObject CrossReturnType; 00559 template<typename OtherDerived> 00560 const CrossReturnType cross(const MatrixBase<OtherDerived>& other) const; 00561 00562 enum { 00563 HNormalized_Size = Direction==Vertical ? internal::traits<ExpressionType>::RowsAtCompileTime 00564 : internal::traits<ExpressionType>::ColsAtCompileTime, 00565 HNormalized_SizeMinusOne = HNormalized_Size==Dynamic ? Dynamic : HNormalized_Size-1 00566 }; 00567 typedef Block<const ExpressionType, 00568 Direction==Vertical ? int(HNormalized_SizeMinusOne) 00569 : int(internal::traits<ExpressionType>::RowsAtCompileTime), 00570 Direction==Horizontal ? int(HNormalized_SizeMinusOne) 00571 : int(internal::traits<ExpressionType>::ColsAtCompileTime)> 00572 HNormalized_Block; 00573 typedef Block<const ExpressionType, 00574 Direction==Vertical ? 1 : int(internal::traits<ExpressionType>::RowsAtCompileTime), 00575 Direction==Horizontal ? 1 : int(internal::traits<ExpressionType>::ColsAtCompileTime)> 00576 HNormalized_Factors; 00577 typedef CwiseBinaryOp<internal::scalar_quotient_op<typename internal::traits<ExpressionType>::Scalar>, 00578 const HNormalized_Block, 00579 const Replicate<HNormalized_Factors, 00580 Direction==Vertical ? HNormalized_SizeMinusOne : 1, 00581 Direction==Horizontal ? HNormalized_SizeMinusOne : 1> > 00582 HNormalizedReturnType; 00583 00584 const HNormalizedReturnType hnormalized() const; 00585 00586 protected: 00587 ExpressionTypeNested m_matrix; 00588 }; 00589 00590 /** \returns a VectorwiseOp wrapper of *this providing additional partial reduction operations 00591 * 00592 * Example: \include MatrixBase_colwise.cpp 00593 * Output: \verbinclude MatrixBase_colwise.out 00594 * 00595 * \sa rowwise(), class VectorwiseOp, \ref TutorialReductionsVisitorsBroadcasting 00596 */ 00597 template<typename Derived> 00598 inline const typename DenseBase<Derived>::ConstColwiseReturnType 00599 DenseBase<Derived>::colwise() const 00600 { 00601 return derived(); 00602 } 00603 00604 /** \returns a writable VectorwiseOp wrapper of *this providing additional partial reduction operations 00605 * 00606 * \sa rowwise(), class VectorwiseOp, \ref TutorialReductionsVisitorsBroadcasting 00607 */ 00608 template<typename Derived> 00609 inline typename DenseBase<Derived>::ColwiseReturnType 00610 DenseBase<Derived>::colwise() 00611 { 00612 return derived(); 00613 } 00614 00615 /** \returns a VectorwiseOp wrapper of *this providing additional partial reduction operations 00616 * 00617 * Example: \include MatrixBase_rowwise.cpp 00618 * Output: \verbinclude MatrixBase_rowwise.out 00619 * 00620 * \sa colwise(), class VectorwiseOp, \ref TutorialReductionsVisitorsBroadcasting 00621 */ 00622 template<typename Derived> 00623 inline const typename DenseBase<Derived>::ConstRowwiseReturnType 00624 DenseBase<Derived>::rowwise() const 00625 { 00626 return derived(); 00627 } 00628 00629 /** \returns a writable VectorwiseOp wrapper of *this providing additional partial reduction operations 00630 * 00631 * \sa colwise(), class VectorwiseOp, \ref TutorialReductionsVisitorsBroadcasting 00632 */ 00633 template<typename Derived> 00634 inline typename DenseBase<Derived>::RowwiseReturnType 00635 DenseBase<Derived>::rowwise() 00636 { 00637 return derived(); 00638 } 00639 00640 } // end namespace Eigen 00641 00642 #endif // EIGEN_PARTIAL_REDUX_H
Generated on Thu Nov 17 2022 22:01:31 by
1.7.2