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