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