Classes | Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Eigen::VectorwiseOp Class Reference

Pseudo expression providing broadcasting and partial reduction operations. More...

#include <ForwardDeclarations.h>

Classes

struct  ExtendedType
 
struct  LpNormReturnType
 
struct  OppositeExtendedType
 
struct  ReduxReturnType
 
struct  ReturnType
 

Public Types

enum  { isVertical = (Direction==Vertical) ? 1 : 0, isHorizontal = (Direction==Horizontal) ? 1 : 0 }
 
enum  { HNormalized_Size, HNormalized_SizeMinusOne = HNormalized_Size==Dynamic ? Dynamic : HNormalized_Size-1 }
 
typedef ReturnType< internal::member_all >::Type AllReturnType
 
typedef ReturnType< internal::member_any >::Type AnyReturnType
 
typedef ReturnType< internal::member_blueNorm, RealScalar >::Type BlueNormReturnType
 
typedef internal::subvector_stl_iterator< const ExpressionType, DirectionType(Direction)> const_iterator
 
typedef internal::subvector_stl_reverse_iterator< const ExpressionType, DirectionType(Direction)> const_reverse_iterator
 
typedef Reverse< const ExpressionType, Direction > ConstReverseReturnType
 
typedef PartialReduxExpr< ExpressionType, internal::member_count< Index, Scalar >, Direction > CountReturnType
 
typedef ExpressionType::PlainObject CrossReturnType
 
typedef internal::ref_selector< ExpressionType >::non_const_type ExpressionTypeNested
 
typedef internal::remove_all< ExpressionTypeNested >::type ExpressionTypeNestedCleaned
 
typedef Block< const ExpressionType, Direction==Vertical ? int(HNormalized_SizeMinusOne) :int(internal::traits< ExpressionType >::RowsAtCompileTime), Direction==Horizontal ? int(HNormalized_SizeMinusOne) :int(internal::traits< ExpressionType >::ColsAtCompileTime)> HNormalized_Block
 
typedef Block< const ExpressionType, Direction==Vertical ? 1 :int(internal::traits< ExpressionType >::RowsAtCompileTime), Direction==Horizontal ? 1 :int(internal::traits< ExpressionType >::ColsAtCompileTime)> HNormalized_Factors
 
typedef CwiseBinaryOp< internal::scalar_quotient_op< typename internal::traits< ExpressionType >::Scalar >, const HNormalized_Block, const Replicate< HNormalized_Factors, Direction==Vertical ? HNormalized_SizeMinusOne :1, Direction==Horizontal ? HNormalized_SizeMinusOne :1 > > HNormalizedReturnType
 
typedef Homogeneous< ExpressionType, Direction > HomogeneousReturnType
 
typedef ReturnType< internal::member_hypotNorm, RealScalar >::Type HypotNormReturnType
 
typedef Eigen::Index Index
 
typedef internal::subvector_stl_iterator< ExpressionType, DirectionType(Direction)> iterator
 
typedef ReturnType< internal::member_maxCoeff >::Type MaxCoeffReturnType
 
typedef ReturnType< internal::member_minCoeff >::Type MinCoeffReturnType
 
typedef CwiseUnaryOp< internal::scalar_sqrt_op< RealScalar >, const SquaredNormReturnTypeNormReturnType
 
typedef ReturnType< internal::member_prod >::Type ProdReturnType
 
typedef ExpressionType::RealScalar RealScalar
 
typedef Replicate< ExpressionType,(isVertical?Dynamic:1),(isHorizontal?Dynamic:1)> ReplicateReturnType
 
typedef internal::subvector_stl_reverse_iterator< ExpressionType, DirectionType(Direction)> reverse_iterator
 
typedef Reverse< ExpressionType, Direction > ReverseReturnType
 
typedef ExpressionType::Scalar Scalar
 
typedef PartialReduxExpr< const CwiseUnaryOp< internal::scalar_abs2_op< Scalar >, const ExpressionTypeNestedCleaned >, internal::member_sum< RealScalar, RealScalar >, Direction > SquaredNormReturnType
 
typedef ReturnType< internal::member_stableNorm, RealScalar >::Type StableNormReturnType
 
typedef ReturnType< internal::member_sum >::Type SumReturnType
 

Public Member Functions

const EIGEN_DEVICE_FUNC ExpressionType & _expression () const
 
const EIGEN_DEVICE_FUNC AllReturnType all () const
 
const EIGEN_DEVICE_FUNC AnyReturnType any () const
 
iterator begin ()
 
const_iterator begin () const
 
const EIGEN_DEVICE_FUNC BlueNormReturnType blueNorm () const
 
const_iterator cbegin () const
 
const_iterator cend () const
 
const EIGEN_DEVICE_FUNC CountReturnType count () const
 
const_reverse_iterator crbegin () const
 
const_reverse_iterator crend () const
 
template<typename OtherDerived >
const EIGEN_DEVICE_FUNC CrossReturnType cross (const MatrixBase< OtherDerived > &other) const
 
typedef EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE (SumReturnType, Scalar, quotient) MeanReturnType
 
iterator end ()
 
const_iterator end () const
 
const EIGEN_DEVICE_FUNC HNormalizedReturnType hnormalized () const
 column or row-wise homogeneous normalization More...
 
EIGEN_DEVICE_FUNC HomogeneousReturnType homogeneous () const
 
const EIGEN_DEVICE_FUNC HypotNormReturnType hypotNorm () const
 
template<int p>
const EIGEN_DEVICE_FUNC LpNormReturnType< p >::Type lpNorm () const
 
const EIGEN_DEVICE_FUNC MaxCoeffReturnType maxCoeff () const
 
const EIGEN_DEVICE_FUNC MeanReturnType mean () const
 
const EIGEN_DEVICE_FUNC MinCoeffReturnType minCoeff () const
 
const EIGEN_DEVICE_FUNC NormReturnType norm () const
 
EIGEN_DEVICE_FUNC void normalize ()
 
EIGEN_DEVICE_FUNC CwiseBinaryOp< internal::scalar_quotient_op< Scalar >, const ExpressionTypeNestedCleaned, const typename OppositeExtendedType< NormReturnType >::Typenormalized () const
 
template<typename OtherDerived >
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC CwiseBinaryOp< internal::scalar_product_op< Scalar >, const ExpressionTypeNestedCleaned, const typename ExtendedType< OtherDerived >::Type > EIGEN_DEVICE_FUNC operator* (const DenseBase< OtherDerived > &other) const
 
template<typename OtherDerived >
EIGEN_DEVICE_FUNC ExpressionType & operator*= (const DenseBase< OtherDerived > &other)
 
template<typename OtherDerived >
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC CwiseBinaryOp< internal::scalar_sum_op< Scalar, typename OtherDerived::Scalar >, const ExpressionTypeNestedCleaned, const typename ExtendedType< OtherDerived >::Typeoperator+ (const DenseBase< OtherDerived > &other) const
 
template<typename OtherDerived >
EIGEN_DEVICE_FUNC ExpressionType & operator+= (const DenseBase< OtherDerived > &other)
 
template<typename OtherDerived >
EIGEN_DEVICE_FUNC CwiseBinaryOp< internal::scalar_difference_op< Scalar, typename OtherDerived::Scalar >, const ExpressionTypeNestedCleaned, const typename ExtendedType< OtherDerived >::Typeoperator- (const DenseBase< OtherDerived > &other) const
 
template<typename OtherDerived >
EIGEN_DEVICE_FUNC ExpressionType & operator-= (const DenseBase< OtherDerived > &other)
 
template<typename OtherDerived >
EIGEN_DEVICE_FUNC CwiseBinaryOp< internal::scalar_quotient_op< Scalar >, const ExpressionTypeNestedCleaned, const typename ExtendedType< OtherDerived >::Typeoperator/ (const DenseBase< OtherDerived > &other) const
 
template<typename OtherDerived >
EIGEN_DEVICE_FUNC ExpressionType & operator/= (const DenseBase< OtherDerived > &other)
 
template<typename OtherDerived >
EIGEN_DEVICE_FUNC ExpressionType & operator= (const DenseBase< OtherDerived > &other)
 
const EIGEN_DEVICE_FUNC ProdReturnType prod () const
 
reverse_iterator rbegin ()
 
const_reverse_iterator rbegin () const
 
template<typename BinaryOp >
const EIGEN_DEVICE_FUNC ReduxReturnType< BinaryOp >::Type redux (const BinaryOp &func=BinaryOp()) const
 
reverse_iterator rend ()
 
const_reverse_iterator rend () const
 
const EIGEN_DEVICE_FUNC ReplicateReturnType replicate (Index factor) const
 
template<int Factor>
const Replicate< ExpressionType, isVertical *Factor+isHorizontal, isHorizontal *Factor+isVertical > EIGEN_DEVICE_FUNC replicate (Index factor=Factor) const
 
EIGEN_DEVICE_FUNC ReverseReturnType reverse ()
 
const EIGEN_DEVICE_FUNC ConstReverseReturnType reverse () const
 
EIGEN_DEVICE_FUNC void reverseInPlace ()
 
const EIGEN_DEVICE_FUNC SquaredNormReturnType squaredNorm () const
 
const EIGEN_DEVICE_FUNC StableNormReturnType stableNorm () const
 
const EIGEN_DEVICE_FUNC SumReturnType sum () const
 
EIGEN_DEVICE_FUNC VectorwiseOp (ExpressionType &matrix)
 

Protected Member Functions

template<typename OtherDerived >
EIGEN_DEVICE_FUNC ExtendedType< OtherDerived >::Type extendedTo (const DenseBase< OtherDerived > &other) const
 
template<typename OtherDerived >
EIGEN_DEVICE_FUNC OppositeExtendedType< OtherDerived >::Type extendedToOpposite (const DenseBase< OtherDerived > &other) const
 
Index redux_length () const
 

Protected Attributes

ExpressionTypeNested m_matrix
 

Detailed Description

Pseudo expression providing broadcasting and partial reduction operations.

Template Parameters
ExpressionTypethe type of the object on which to do partial reductions
Directionindicates whether to operate on columns (Vertical) or rows (Horizontal)

This class represents a pseudo expression with broadcasting and partial reduction features. It is the return type of DenseBase::colwise() and DenseBase::rowwise() and most of the time this is the only way it is explicitly used.

To understand the logic of rowwise/colwise expression, let's consider a generic case A.colwise().foo() where foo is any method of VectorwiseOp. This expression is equivalent to applying foo() to each column of A and then re-assemble the outputs in a matrix expression:

[A.col(0).foo(), A.col(1).foo(), ..., A.col(A.cols()-1).foo()]

Example:

Matrix3d m = Matrix3d::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the sum of each column:" << endl << m.colwise().sum() << endl;
cout << "Here is the maximum absolute value of each column:"
<< endl << m.cwiseAbs().colwise().maxCoeff() << endl;

Output:

The begin() and end() methods are obviously exceptions to the previous rule as they return STL-compatible begin/end iterators to the rows or columns of the nested expression. Typical use cases include for-range-loop and calls to STL algorithms:

Example:

Matrix3i m = Matrix3i::Random();
cout << "Here is the initial matrix m:" << endl << m << endl;
int i = -1;
for(auto c: m.colwise()) {
c *= i;
++i;
}
cout << "Here is the matrix m after the for-range-loop:" << endl << m << endl;
auto cols = m.colwise();
auto it = std::find_if(cols.cbegin(), cols.cend(),
[](Matrix3i::ConstColXpr x) { return x.squaredNorm() == 0; });
cout << "The first empty column is: " << distance(cols.cbegin(),it) << endl;

Output:

For a partial reduction on an empty input, some rules apply. For the sake of clarity, let's consider a vertical reduction:

See also
DenseBase::colwise(), DenseBase::rowwise(), class PartialReduxExpr

Definition at line 264 of file ForwardDeclarations.h.

Member Typedef Documentation

◆ AllReturnType

typedef ReturnType<internal::member_all>::Type Eigen::VectorwiseOp::AllReturnType

Definition at line 352 of file VectorwiseOp.h.

◆ AnyReturnType

typedef ReturnType<internal::member_any>::Type Eigen::VectorwiseOp::AnyReturnType

Definition at line 353 of file VectorwiseOp.h.

◆ BlueNormReturnType

Definition at line 347 of file VectorwiseOp.h.

◆ const_iterator

Definition at line 283 of file VectorwiseOp.h.

◆ const_reverse_iterator

Definition at line 285 of file VectorwiseOp.h.

◆ ConstReverseReturnType

typedef Reverse<const ExpressionType, Direction> Eigen::VectorwiseOp::ConstReverseReturnType

Definition at line 356 of file VectorwiseOp.h.

◆ CountReturnType

typedef PartialReduxExpr<ExpressionType, internal::member_count<Index,Scalar>, Direction> Eigen::VectorwiseOp::CountReturnType

Definition at line 354 of file VectorwiseOp.h.

◆ CrossReturnType

typedef ExpressionType::PlainObject Eigen::VectorwiseOp::CrossReturnType

Definition at line 712 of file VectorwiseOp.h.

◆ ExpressionTypeNested

typedef internal::ref_selector<ExpressionType>::non_const_type Eigen::VectorwiseOp::ExpressionTypeNested

Definition at line 193 of file VectorwiseOp.h.

◆ ExpressionTypeNestedCleaned

Definition at line 194 of file VectorwiseOp.h.

◆ HNormalized_Block

typedef Block<const ExpressionType, Direction==Vertical ? int(HNormalized_SizeMinusOne) : int(internal::traits<ExpressionType>::RowsAtCompileTime), Direction==Horizontal ? int(HNormalized_SizeMinusOne) : int(internal::traits<ExpressionType>::ColsAtCompileTime)> Eigen::VectorwiseOp::HNormalized_Block

Definition at line 727 of file VectorwiseOp.h.

◆ HNormalized_Factors

typedef Block<const ExpressionType, Direction==Vertical ? 1 : int(internal::traits<ExpressionType>::RowsAtCompileTime), Direction==Horizontal ? 1 : int(internal::traits<ExpressionType>::ColsAtCompileTime)> Eigen::VectorwiseOp::HNormalized_Factors

Definition at line 731 of file VectorwiseOp.h.

◆ HNormalizedReturnType

Definition at line 737 of file VectorwiseOp.h.

◆ HomogeneousReturnType

typedef Homogeneous<ExpressionType,Direction> Eigen::VectorwiseOp::HomogeneousReturnType

Definition at line 708 of file VectorwiseOp.h.

◆ HypotNormReturnType

Definition at line 349 of file VectorwiseOp.h.

◆ Index

Deprecated:
since Eigen 3.3

Definition at line 192 of file VectorwiseOp.h.

◆ iterator

Definition at line 282 of file VectorwiseOp.h.

◆ MaxCoeffReturnType

typedef ReturnType<internal::member_maxCoeff>::Type Eigen::VectorwiseOp::MaxCoeffReturnType

Definition at line 344 of file VectorwiseOp.h.

◆ MinCoeffReturnType

typedef ReturnType<internal::member_minCoeff>::Type Eigen::VectorwiseOp::MinCoeffReturnType

Definition at line 343 of file VectorwiseOp.h.

◆ NormReturnType

Definition at line 346 of file VectorwiseOp.h.

◆ ProdReturnType

typedef ReturnType<internal::member_prod>::Type Eigen::VectorwiseOp::ProdReturnType

Definition at line 355 of file VectorwiseOp.h.

◆ RealScalar

typedef ExpressionType::RealScalar Eigen::VectorwiseOp::RealScalar

Definition at line 191 of file VectorwiseOp.h.

◆ ReplicateReturnType

typedef Replicate<ExpressionType,(isVertical?Dynamic:1),(isHorizontal?Dynamic:1)> Eigen::VectorwiseOp::ReplicateReturnType

Definition at line 552 of file VectorwiseOp.h.

◆ reverse_iterator

Definition at line 284 of file VectorwiseOp.h.

◆ ReverseReturnType

typedef Reverse<ExpressionType, Direction> Eigen::VectorwiseOp::ReverseReturnType

Definition at line 357 of file VectorwiseOp.h.

◆ Scalar

typedef ExpressionType::Scalar Eigen::VectorwiseOp::Scalar

Definition at line 190 of file VectorwiseOp.h.

◆ SquaredNormReturnType

Definition at line 345 of file VectorwiseOp.h.

◆ StableNormReturnType

Definition at line 348 of file VectorwiseOp.h.

◆ SumReturnType

typedef ReturnType<internal::member_sum>::Type Eigen::VectorwiseOp::SumReturnType

Definition at line 350 of file VectorwiseOp.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
isVertical 
isHorizontal 

Definition at line 213 of file VectorwiseOp.h.

◆ anonymous enum

anonymous enum
Enumerator
HNormalized_Size 
HNormalized_SizeMinusOne 

Definition at line 717 of file VectorwiseOp.h.

Constructor & Destructor Documentation

◆ VectorwiseOp()

EIGEN_DEVICE_FUNC Eigen::VectorwiseOp::VectorwiseOp ( ExpressionType &  matrix)
inlineexplicit

Definition at line 268 of file VectorwiseOp.h.

Member Function Documentation

◆ _expression()

const EIGEN_DEVICE_FUNC ExpressionType& Eigen::VectorwiseOp::_expression ( ) const
inline

Definition at line 272 of file VectorwiseOp.h.

◆ all()

const EIGEN_DEVICE_FUNC AllReturnType Eigen::VectorwiseOp::all ( ) const
inline
Returns
a row (or column) vector expression representing whether all coefficients of each respective column (or row) are true. This expression can be assigned to a vector with entries of type bool.
See also
DenseBase::all()

Definition at line 496 of file VectorwiseOp.h.

◆ any()

const EIGEN_DEVICE_FUNC AnyReturnType Eigen::VectorwiseOp::any ( ) const
inline
Returns
a row (or column) vector expression representing whether at least one coefficient of each respective column (or row) is true. This expression can be assigned to a vector with entries of type bool.
See also
DenseBase::any()

Definition at line 505 of file VectorwiseOp.h.

◆ begin() [1/2]

iterator Eigen::VectorwiseOp::begin ( )
inline

returns an iterator to the first row (rowwise) or column (colwise) of the nested expression.

See also
end(), cbegin()

Definition at line 291 of file VectorwiseOp.h.

◆ begin() [2/2]

const_iterator Eigen::VectorwiseOp::begin ( ) const
inline

const version of begin()

Definition at line 293 of file VectorwiseOp.h.

◆ blueNorm()

const EIGEN_DEVICE_FUNC BlueNormReturnType Eigen::VectorwiseOp::blueNorm ( ) const
inline
Returns
a row (or column) vector expression of the norm of each column (or row) of the referenced expression, using Blue's algorithm. This is a vector with real entries, even if the original matrix has complex entries.
See also
DenseBase::blueNorm()

Definition at line 446 of file VectorwiseOp.h.

◆ cbegin()

const_iterator Eigen::VectorwiseOp::cbegin ( ) const
inline

const version of begin()

Definition at line 295 of file VectorwiseOp.h.

◆ cend()

const_iterator Eigen::VectorwiseOp::cend ( ) const
inline

const version of end()

Definition at line 313 of file VectorwiseOp.h.

◆ count()

const EIGEN_DEVICE_FUNC CountReturnType Eigen::VectorwiseOp::count ( ) const
inline
Returns
a row (or column) vector expression representing the number of true coefficients of each respective column (or row). This expression can be assigned to a vector whose entries have the same type as is used to index entries of the original matrix; for dense matrices, this is std::ptrdiff_t .

Example:

Matrix3d m = Matrix3d::Random();
cout << "Here is the matrix m:" << endl << m << endl;
Matrix<ptrdiff_t, 3, 1> res = (m.array() >= 0.5).rowwise().count();
cout << "Here is the count of elements larger or equal than 0.5 of each row:" << endl;
cout << res << endl;

Output:

See also
DenseBase::count()

Definition at line 518 of file VectorwiseOp.h.

◆ crbegin()

const_reverse_iterator Eigen::VectorwiseOp::crbegin ( ) const
inline

const version of rbegin()

Definition at line 304 of file VectorwiseOp.h.

◆ crend()

const_reverse_iterator Eigen::VectorwiseOp::crend ( ) const
inline

const version of rend()

Definition at line 322 of file VectorwiseOp.h.

◆ EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE()

typedef Eigen::VectorwiseOp::EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE ( SumReturnType  ,
Scalar  ,
quotient   
)

◆ end() [1/2]

iterator Eigen::VectorwiseOp::end ( )
inline

returns an iterator to the row (resp. column) following the last row (resp. column) of the nested expression

See also
begin(), cend()

Definition at line 309 of file VectorwiseOp.h.

◆ end() [2/2]

const_iterator Eigen::VectorwiseOp::end ( ) const
inline

const version of end()

Definition at line 311 of file VectorwiseOp.h.

◆ extendedTo()

template<typename OtherDerived >
EIGEN_DEVICE_FUNC ExtendedType<OtherDerived>::Type Eigen::VectorwiseOp::extendedTo ( const DenseBase< OtherDerived > &  other) const
inlineprotected

Definition at line 231 of file VectorwiseOp.h.

◆ extendedToOpposite()

template<typename OtherDerived >
EIGEN_DEVICE_FUNC OppositeExtendedType<OtherDerived>::Type Eigen::VectorwiseOp::extendedToOpposite ( const DenseBase< OtherDerived > &  other) const
inlineprotected

Definition at line 254 of file VectorwiseOp.h.

◆ hypotNorm()

const EIGEN_DEVICE_FUNC HypotNormReturnType Eigen::VectorwiseOp::hypotNorm ( ) const
inline
Returns
a row (or column) vector expression of the norm of each column (or row) of the referenced expression, avoiding underflow and overflow using a concatenation of hypot() calls. This is a vector with real entries, even if the original matrix has complex entries.
See also
DenseBase::hypotNorm()

Definition at line 468 of file VectorwiseOp.h.

◆ lpNorm()

template<int p>
const EIGEN_DEVICE_FUNC LpNormReturnType<p>::Type Eigen::VectorwiseOp::lpNorm ( ) const
inline
Returns
a row (or column) vector expression of the norm of each column (or row) of the referenced expression. This is a vector with real entries, even if the original matrix has complex entries.

Example:

Matrix3d m = Matrix3d::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the norm of each column:" << endl << m.colwise().norm() << endl;

Output:

See also
DenseBase::norm()

Definition at line 435 of file VectorwiseOp.h.

◆ maxCoeff()

const EIGEN_DEVICE_FUNC MaxCoeffReturnType Eigen::VectorwiseOp::maxCoeff ( ) const
inline
Returns
a row (or column) vector expression of the largest coefficient of each column (or row) of the referenced expression.
Warning
the size along the reduction direction must be strictly positive, otherwise an assertion is triggered.
the result is undefined if *this contains NaN.

Example:

Matrix3d m = Matrix3d::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the maximum of each column:" << endl << m.colwise().maxCoeff() << endl;

Output:

See also
DenseBase::maxCoeff()

Definition at line 395 of file VectorwiseOp.h.

◆ mean()

const EIGEN_DEVICE_FUNC MeanReturnType Eigen::VectorwiseOp::mean ( ) const
inline
Returns
a row (or column) vector expression of the mean of each column (or row) of the referenced expression.
See also
DenseBase::mean()

Definition at line 487 of file VectorwiseOp.h.

◆ minCoeff()

const EIGEN_DEVICE_FUNC MinCoeffReturnType Eigen::VectorwiseOp::minCoeff ( ) const
inline
Returns
a row (or column) vector expression of the smallest coefficient of each column (or row) of the referenced expression.
Warning
the size along the reduction direction must be strictly positive, otherwise an assertion is triggered.
the result is undefined if *this contains NaN.

Example:

Matrix3d m = Matrix3d::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the minimum of each column:" << endl << m.colwise().minCoeff() << endl;

Output:

See also
DenseBase::minCoeff()

Definition at line 376 of file VectorwiseOp.h.

◆ norm()

const EIGEN_DEVICE_FUNC NormReturnType Eigen::VectorwiseOp::norm ( ) const
inline
Returns
a row (or column) vector expression of the norm of each column (or row) of the referenced expression. This is a vector with real entries, even if the original matrix has complex entries.

Example:

Matrix3d m = Matrix3d::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the norm of each column:" << endl << m.colwise().norm() << endl;

Output:

See also
DenseBase::norm()

Definition at line 422 of file VectorwiseOp.h.

◆ normalize()

EIGEN_DEVICE_FUNC void Eigen::VectorwiseOp::normalize ( )
inline

Normalize in-place each row or columns of the referenced matrix.

See also
MatrixBase::normalize(), normalized()

Definition at line 700 of file VectorwiseOp.h.

◆ normalized()

Returns
an expression where each column (or row) of the referenced matrix are normalized. The referenced matrix is not modified.
See also
MatrixBase::normalized(), normalize()

Definition at line 694 of file VectorwiseOp.h.

◆ operator*()

template<typename OtherDerived >
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC CwiseBinaryOp<internal::scalar_product_op<Scalar>, const ExpressionTypeNestedCleaned, const typename ExtendedType<OtherDerived>::Type> EIGEN_DEVICE_FUNC Eigen::VectorwiseOp::operator* ( const DenseBase< OtherDerived > &  other) const
inline

Returns the expression where each subvector is the product of the vector other by the corresponding subvector of *this

Definition at line 663 of file VectorwiseOp.h.

◆ operator*=()

template<typename OtherDerived >
EIGEN_DEVICE_FUNC ExpressionType& Eigen::VectorwiseOp::operator*= ( const DenseBase< OtherDerived > &  other)
inline

Multiples each subvector of *this by the vector other

Definition at line 610 of file VectorwiseOp.h.

◆ operator+()

template<typename OtherDerived >
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC CwiseBinaryOp<internal::scalar_sum_op<Scalar,typename OtherDerived::Scalar>, const ExpressionTypeNestedCleaned, const typename ExtendedType<OtherDerived>::Type> Eigen::VectorwiseOp::operator+ ( const DenseBase< OtherDerived > &  other) const
inline

Returns the expression of the sum of the vector other to each subvector of *this

Definition at line 636 of file VectorwiseOp.h.

◆ operator+=()

template<typename OtherDerived >
EIGEN_DEVICE_FUNC ExpressionType& Eigen::VectorwiseOp::operator+= ( const DenseBase< OtherDerived > &  other)
inline

Adds the vector other to each subvector of *this

Definition at line 590 of file VectorwiseOp.h.

◆ operator-()

template<typename OtherDerived >
EIGEN_DEVICE_FUNC CwiseBinaryOp<internal::scalar_difference_op<Scalar,typename OtherDerived::Scalar>, const ExpressionTypeNestedCleaned, const typename ExtendedType<OtherDerived>::Type> Eigen::VectorwiseOp::operator- ( const DenseBase< OtherDerived > &  other) const
inline

Returns the expression of the difference between each subvector of *this and the vector other

Definition at line 649 of file VectorwiseOp.h.

◆ operator-=()

template<typename OtherDerived >
EIGEN_DEVICE_FUNC ExpressionType& Eigen::VectorwiseOp::operator-= ( const DenseBase< OtherDerived > &  other)
inline

Substracts the vector other to each subvector of *this

Definition at line 600 of file VectorwiseOp.h.

◆ operator/()

template<typename OtherDerived >
EIGEN_DEVICE_FUNC CwiseBinaryOp<internal::scalar_quotient_op<Scalar>, const ExpressionTypeNestedCleaned, const typename ExtendedType<OtherDerived>::Type> Eigen::VectorwiseOp::operator/ ( const DenseBase< OtherDerived > &  other) const
inline

Returns the expression where each subvector is the quotient of the corresponding subvector of *this by the vector other

Definition at line 678 of file VectorwiseOp.h.

◆ operator/=()

template<typename OtherDerived >
EIGEN_DEVICE_FUNC ExpressionType& Eigen::VectorwiseOp::operator/= ( const DenseBase< OtherDerived > &  other)
inline

Divides each subvector of *this by the vector other

Definition at line 622 of file VectorwiseOp.h.

◆ operator=()

template<typename OtherDerived >
EIGEN_DEVICE_FUNC ExpressionType& Eigen::VectorwiseOp::operator= ( const DenseBase< OtherDerived > &  other)
inline

Copies the vector other to each subvector of *this

Definition at line 579 of file VectorwiseOp.h.

◆ prod()

const EIGEN_DEVICE_FUNC ProdReturnType Eigen::VectorwiseOp::prod ( ) const
inline
Returns
a row (or column) vector expression of the product of each column (or row) of the referenced expression.

Example:

Matrix3d m = Matrix3d::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the product of each row:" << endl << m.rowwise().prod() << endl;

Output:

See also
DenseBase::prod()

Definition at line 529 of file VectorwiseOp.h.

◆ rbegin() [1/2]

reverse_iterator Eigen::VectorwiseOp::rbegin ( )
inline

returns a reverse iterator to the last row (rowwise) or column (colwise) of the nested expression.

See also
rend(), crbegin()

Definition at line 300 of file VectorwiseOp.h.

◆ rbegin() [2/2]

const_reverse_iterator Eigen::VectorwiseOp::rbegin ( ) const
inline

const version of rbegin()

Definition at line 302 of file VectorwiseOp.h.

◆ redux()

template<typename BinaryOp >
const EIGEN_DEVICE_FUNC ReduxReturnType<BinaryOp>::Type Eigen::VectorwiseOp::redux ( const BinaryOp &  func = BinaryOp()) const
inline
Returns
a row or column vector expression of *this reduxed by func

The template parameter BinaryOp is the type of the functor of the custom redux operator. Note that func must be an associative operator.

Warning
the size along the reduction direction must be strictly positive, otherwise an assertion is triggered.
See also
class VectorwiseOp, DenseBase::colwise(), DenseBase::rowwise()

Definition at line 337 of file VectorwiseOp.h.

◆ redux_length()

Index Eigen::VectorwiseOp::redux_length ( ) const
inlineprotected

Definition at line 747 of file VectorwiseOp.h.

◆ rend() [1/2]

reverse_iterator Eigen::VectorwiseOp::rend ( )
inline

returns a reverse iterator to the row (resp. column) before the first row (resp. column) of the nested expression

See also
begin(), cend()

Definition at line 318 of file VectorwiseOp.h.

◆ rend() [2/2]

const_reverse_iterator Eigen::VectorwiseOp::rend ( ) const
inline

const version of rend()

Definition at line 320 of file VectorwiseOp.h.

◆ replicate() [1/2]

const EIGEN_DEVICE_FUNC VectorwiseOp< ExpressionType, Direction >::ReplicateReturnType Eigen::VectorwiseOp::replicate ( Index  factor) const
Returns
an expression of the replication of each column (or row) of *this

Example:

Vector3i v = Vector3i::Random();
cout << "Here is the vector v:" << endl << v << endl;
cout << "v.rowwise().replicate(5) = ..." << endl;
cout << v.rowwise().replicate(5) << endl;

Output:

See also
VectorwiseOp::replicate(), DenseBase::replicate(), class Replicate

Definition at line 134 of file Replicate.h.

◆ replicate() [2/2]

template<int Factor>
const Replicate<ExpressionType,isVertical*Factor+isHorizontal,isHorizontal*Factor+isVertical> EIGEN_DEVICE_FUNC Eigen::VectorwiseOp::replicate ( Index  factor = Factor) const
inline
Returns
an expression of the replication of each column (or row) of *this

Example:

MatrixXi m = MatrixXi::Random(2,3);
cout << "Here is the matrix m:" << endl << m << endl;
cout << "m.colwise().replicate<3>() = ..." << endl;
cout << m.colwise().replicate<3>() << endl;

Output:

See also
VectorwiseOp::replicate(Index), DenseBase::replicate(), class Replicate

Definition at line 568 of file VectorwiseOp.h.

◆ reverse() [1/2]

EIGEN_DEVICE_FUNC ReverseReturnType Eigen::VectorwiseOp::reverse ( )
inline
Returns
a writable matrix expression where each column (or row) are reversed.
See also
reverse() const

Definition at line 549 of file VectorwiseOp.h.

◆ reverse() [2/2]

const EIGEN_DEVICE_FUNC ConstReverseReturnType Eigen::VectorwiseOp::reverse ( ) const
inline
Returns
a matrix expression where each column (or row) are reversed.

Example:

MatrixXi m = MatrixXi::Random(3,4);
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the rowwise reverse of m:" << endl << m.rowwise().reverse() << endl;
cout << "Here is the colwise reverse of m:" << endl << m.colwise().reverse() << endl;
cout << "Here is the coefficient (1,0) in the rowise reverse of m:" << endl
<< m.rowwise().reverse()(1,0) << endl;
cout << "Let us overwrite this coefficient with the value 4." << endl;
//m.colwise().reverse()(1,0) = 4;
cout << "Now the matrix m is:" << endl << m << endl;

Output:

See also
DenseBase::reverse()

Definition at line 541 of file VectorwiseOp.h.

◆ reverseInPlace()

EIGEN_DEVICE_FUNC void Eigen::VectorwiseOp::reverseInPlace ( )
inline

This is the "in place" version of VectorwiseOp::reverse: it reverses each column or row of *this.

In most cases it is probably better to simply use the reversed expression of a matrix. However, when reversing the matrix data itself is really needed, then this "in-place" version is probably the right choice because it provides the following additional benefits:

  • less error prone: doing the same operation with .reverse() requires special care:
    m = m.reverse().eval();
  • this API enables reverse operations without the need for a temporary
See also
DenseBase::reverseInPlace(), reverse()

Definition at line 210 of file Reverse.h.

◆ squaredNorm()

const EIGEN_DEVICE_FUNC SquaredNormReturnType Eigen::VectorwiseOp::squaredNorm ( ) const
inline
Returns
a row (or column) vector expression of the squared norm of each column (or row) of the referenced expression. This is a vector with real entries, even if the original matrix has complex entries.

Example:

Matrix3d m = Matrix3d::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the square norm of each row:" << endl << m.rowwise().squaredNorm() << endl;

Output:

See also
DenseBase::squaredNorm()

Definition at line 410 of file VectorwiseOp.h.

◆ stableNorm()

const EIGEN_DEVICE_FUNC StableNormReturnType Eigen::VectorwiseOp::stableNorm ( ) const
inline
Returns
a row (or column) vector expression of the norm of each column (or row) of the referenced expression, avoiding underflow and overflow. This is a vector with real entries, even if the original matrix has complex entries.
See also
DenseBase::stableNorm()

Definition at line 457 of file VectorwiseOp.h.

◆ sum()

const EIGEN_DEVICE_FUNC SumReturnType Eigen::VectorwiseOp::sum ( ) const
inline
Returns
a row (or column) vector expression of the sum of each column (or row) of the referenced expression.

Example:

Matrix3d m = Matrix3d::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the sum of each row:" << endl << m.rowwise().sum() << endl;

Output:

See also
DenseBase::sum()

Definition at line 479 of file VectorwiseOp.h.

Member Data Documentation

◆ m_matrix

ExpressionTypeNested Eigen::VectorwiseOp::m_matrix
protected

Definition at line 751 of file VectorwiseOp.h.


The documentation for this class was generated from the following files:
c
Scalar Scalar * c
Definition: benchVecAdd.cpp:17
x
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition: gnuplot_common_settings.hh:12
res
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition: PartialRedux_count.cpp:3
A
Definition: test_numpy_dtypes.cpp:298
m
Matrix3f m
Definition: AngleAxis_mimic_euler.cpp:1
ConstColXpr
const typedef Block< const Derived, internal::traits< Derived >::RowsAtCompileTime, 1, !IsRowMajor > ConstColXpr
Definition: BlockMethods.h:15
v
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
gtsam::distance
Double_ distance(const OrientedPlane3_ &p)
Definition: slam/expressions.h:117
cols
int cols
Definition: Tutorial_commainit_02.cpp:1
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9


gtsam
Author(s):
autogenerated on Thu Dec 19 2024 04:11:43