* rename PartialRedux to VectorwiseOp

* add VectorwiseOp's +, -, +=, -= operators
This commit is contained in:
Gael Guennebaud 2009-06-10 11:20:30 +02:00
parent f3fd7fd22b
commit 627595ad19
13 changed files with 174 additions and 52 deletions

View File

@ -32,7 +32,7 @@ namespace Eigen {
#include "src/Array/Functors.h"
#include "src/Array/BooleanRedux.h"
#include "src/Array/Select.h"
#include "src/Array/PartialRedux.h"
#include "src/Array/VectorwiseOp.h"
#include "src/Array/Random.h"
#include "src/Array/Norms.h"
#include "src/Array/Replicate.h"

View File

@ -99,7 +99,7 @@ template<typename MatrixType,int RowFactor,int ColFactor> class Replicate
* Example: \include MatrixBase_replicate.cpp
* Output: \verbinclude MatrixBase_replicate.out
*
* \sa PartialRedux::replicate(), MatrixBase::replicate(int,int), class Replicate
* \sa VectorwiseOp::replicate(), MatrixBase::replicate(int,int), class Replicate
*/
template<typename Derived>
template<int RowFactor, int ColFactor>
@ -115,7 +115,7 @@ MatrixBase<Derived>::replicate() const
* Example: \include MatrixBase_replicate_int_int.cpp
* Output: \verbinclude MatrixBase_replicate_int_int.out
*
* \sa PartialRedux::replicate(), MatrixBase::replicate<int,int>(), class Replicate
* \sa VectorwiseOp::replicate(), MatrixBase::replicate<int,int>(), class Replicate
*/
template<typename Derived>
inline const Replicate<Derived,Dynamic,Dynamic>
@ -130,11 +130,11 @@ MatrixBase<Derived>::replicate(int rowFactor,int colFactor) const
* Example: \include DirectionWise_replicate_int.cpp
* Output: \verbinclude DirectionWise_replicate_int.out
*
* \sa PartialRedux::replicate(), MatrixBase::replicate(), class Replicate
* \sa VectorwiseOp::replicate(), MatrixBase::replicate(), class Replicate
*/
template<typename ExpressionType, int Direction>
const Replicate<ExpressionType,(Direction==Vertical?Dynamic:1),(Direction==Horizontal?Dynamic:1)>
PartialRedux<ExpressionType,Direction>::replicate(int factor) const
VectorwiseOp<ExpressionType,Direction>::replicate(int factor) const
{
return Replicate<ExpressionType,Direction==Vertical?Dynamic:1,Direction==Horizontal?Dynamic:1>
(_expression(),Direction==Vertical?factor:1,Direction==Horizontal?factor:1);
@ -146,12 +146,12 @@ PartialRedux<ExpressionType,Direction>::replicate(int factor) const
* Example: \include DirectionWise_replicate.cpp
* Output: \verbinclude DirectionWise_replicate.out
*
* \sa PartialRedux::replicate(int), MatrixBase::replicate(), class Replicate
* \sa VectorwiseOp::replicate(int), MatrixBase::replicate(), class Replicate
*/
template<typename ExpressionType, int Direction>
template<int Factor>
const Replicate<ExpressionType,(Direction==Vertical?Factor:1),(Direction==Horizontal?Factor:1)>
PartialRedux<ExpressionType,Direction>::replicate(int factor) const
VectorwiseOp<ExpressionType,Direction>::replicate(int factor) const
{
return Replicate<ExpressionType,Direction==Vertical?Factor:1,Direction==Horizontal?Factor:1>
(_expression(),Direction==Vertical?factor:1,Direction==Horizontal?factor:1);

View File

@ -36,10 +36,10 @@
* \param MatrixType the type of the object of which we are taking the reverse
*
* This class represents an expression of the reverse of a vector.
* It is the return type of MatrixBase::reverse() and PartialRedux::reverse()
* It is the return type of MatrixBase::reverse() and VectorwiseOp::reverse()
* and most of the time this is the only way it is used.
*
* \sa MatrixBase::reverse(), PartialRedux::reverse()
* \sa MatrixBase::reverse(), VectorwiseOp::reverse()
*/
template<typename MatrixType, int Direction>
struct ei_traits<Reverse<MatrixType, Direction> >

View File

@ -37,10 +37,10 @@
* \param Direction indicates the direction of the redux (Vertical or Horizontal)
*
* This class represents an expression of a partial redux operator of a matrix.
* It is the return type of PartialRedux functions,
* It is the return type of some VectorwiseOp functions,
* and most of the time this is the only way it is used.
*
* \sa class PartialRedux
* \sa class VectorwiseOp
*/
template< typename MatrixType, typename MemberOp, int Direction>
@ -139,7 +139,7 @@ struct ei_member_redux {
/** \array_module \ingroup Array_Module
*
* \class PartialRedux
* \class VectorwiseOp
*
* \brief Pseudo expression providing partial reduction operations
*
@ -155,7 +155,7 @@ struct ei_member_redux {
*
* \sa MatrixBase::colwise(), MatrixBase::rowwise(), class PartialReduxExpr
*/
template<typename ExpressionType, int Direction> class PartialRedux
template<typename ExpressionType, int Direction> class VectorwiseOp
{
public:
@ -179,7 +179,45 @@ template<typename ExpressionType, int Direction> class PartialRedux
> Type;
};
inline PartialRedux(const ExpressionType& matrix) : m_matrix(matrix) {}
protected:
/** \internal
* \returns the i-th subvector according to the \c Direction */
typedef typename ei_meta_if<Direction==Vertical,
typename ExpressionType::ColXpr,
typename ExpressionType::RowXpr>::ret SubVector;
SubVector subVector(int i)
{
return SubVector(m_matrix.derived(),i);
}
/** \internal
* \returns the number of subvectors in the direction \c Direction */
int subVectors() const
{ return Direction==Vertical?m_matrix.cols():m_matrix.rows(); }
template<typename OtherDerived> struct ExtendedType {
typedef Replicate<OtherDerived,
Direction==Vertical ? 1 : ExpressionType::RowsAtCompileTime,
Direction==Horizontal ? 1 : ExpressionType::ColsAtCompileTime> Type;
};
/** \internal
* Replicates a vector to match the size of \c *this */
template<typename OtherDerived>
typename ExtendedType<OtherDerived>::Type
extendedTo(const MatrixBase<OtherDerived>& other) const
{
EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived);
return typename ExtendedType<OtherDerived>::Type
(other.derived(),
Direction==Vertical ? 1 : m_matrix.rows(),
Direction==Horizontal ? 1 : m_matrix.cols());
}
public:
inline VectorwiseOp(const ExpressionType& matrix) : m_matrix(matrix) {}
/** \internal */
inline const ExpressionType& _expression() const { return m_matrix; }
@ -292,6 +330,48 @@ template<typename ExpressionType, int Direction> class PartialRedux
const Replicate<ExpressionType,(Direction==Vertical?Factor:1),(Direction==Horizontal?Factor:1)>
replicate(int factor = Factor) const;
/////////// Artithmetic operators ///////////
/** Adds the vector \a other to each subvector of \c *this */
template<typename OtherDerived>
ExpressionType& operator+=(const MatrixBase<OtherDerived>& other)
{
for(int j=0; j<subVectors(); ++j)
subVector(j) += other;
return const_cast<ExpressionType&>(m_matrix);
}
/** Substracts the vector \a other to each subvector of \c *this */
template<typename OtherDerived>
ExpressionType& operator-=(const MatrixBase<OtherDerived>& other)
{
for(int j=0; j<subVectors(); ++j)
subVector(j) -= other;
return const_cast<ExpressionType&>(m_matrix);
}
/** Returns the expression of the sum of the vector \a other to each subvector of \c *this */
template<typename OtherDerived>
CwiseBinaryOp<ei_scalar_sum_op<Scalar>,
ExpressionType,
NestByValue<typename ExtendedType<OtherDerived>::Type> >
operator+(const MatrixBase<OtherDerived>& other) const
{
EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived);
return m_matrix + extendedTo(other).nestByValue();
}
/** Returns the expression of the difference between each subvector of \c *this and the vector \a other */
template<typename OtherDerived>
CwiseBinaryOp<ei_scalar_difference_op<Scalar>,
ExpressionType,
NestByValue<typename ExtendedType<OtherDerived>::Type> >
operator-(const MatrixBase<OtherDerived>& other) const
{
EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived);
return m_matrix - extendedTo(other).nestByValue();
}
/////////// Geometry module ///////////
const Homogeneous<ExpressionType,Direction> homogeneous() const;
@ -330,15 +410,15 @@ template<typename ExpressionType, int Direction> class PartialRedux
/** \array_module
*
* \returns a PartialRedux wrapper of *this providing additional partial reduction operations
* \returns a VectorwiseOp wrapper of *this providing additional partial reduction operations
*
* Example: \include MatrixBase_colwise.cpp
* Output: \verbinclude MatrixBase_colwise.out
*
* \sa rowwise(), class PartialRedux
* \sa rowwise(), class VectorwiseOp
*/
template<typename Derived>
inline const PartialRedux<Derived,Vertical>
inline const VectorwiseOp<Derived,Vertical>
MatrixBase<Derived>::colwise() const
{
return derived();
@ -346,31 +426,57 @@ MatrixBase<Derived>::colwise() const
/** \array_module
*
* \returns a PartialRedux wrapper of *this providing additional partial reduction operations
* \returns a writable VectorwiseOp wrapper of *this providing additional partial reduction operations
*
* \sa rowwise(), class VectorwiseOp
*/
template<typename Derived>
inline VectorwiseOp<Derived,Vertical>
MatrixBase<Derived>::colwise()
{
return derived();
}
/** \array_module
*
* \returns a VectorwiseOp wrapper of *this providing additional partial reduction operations
*
* Example: \include MatrixBase_rowwise.cpp
* Output: \verbinclude MatrixBase_rowwise.out
*
* \sa colwise(), class PartialRedux
* \sa colwise(), class VectorwiseOp
*/
template<typename Derived>
inline const PartialRedux<Derived,Horizontal>
inline const VectorwiseOp<Derived,Horizontal>
MatrixBase<Derived>::rowwise() const
{
return derived();
}
/** \array_module
*
* \returns a writable VectorwiseOp wrapper of *this providing additional partial reduction operations
*
* \sa colwise(), class VectorwiseOp
*/
template<typename Derived>
inline VectorwiseOp<Derived,Horizontal>
MatrixBase<Derived>::rowwise()
{
return derived();
}
/** \returns a row or column vector expression of \c *this reduxed by \a func
*
* The template parameter \a BinaryOp is the type of the functor
* of the custom redux operator. Note that func must be an associative operator.
*
* \sa class PartialRedux, MatrixBase::colwise(), MatrixBase::rowwise()
* \sa class VectorwiseOp, MatrixBase::colwise(), MatrixBase::rowwise()
*/
template<typename ExpressionType, int Direction>
template<typename BinaryOp>
const typename PartialRedux<ExpressionType,Direction>::template ReduxReturnType<BinaryOp>::Type
PartialRedux<ExpressionType,Direction>::redux(const BinaryOp& func) const
const typename VectorwiseOp<ExpressionType,Direction>::template ReduxReturnType<BinaryOp>::Type
VectorwiseOp<ExpressionType,Direction>::redux(const BinaryOp& func) const
{
return typename ReduxReturnType<BinaryOp>::Type(_expression(), func);
}

View File

@ -30,7 +30,7 @@
/** \internal
* \brief Template functor to compute the sum of two scalars
*
* \sa class CwiseBinaryOp, MatrixBase::operator+, class PartialRedux, MatrixBase::sum()
* \sa class CwiseBinaryOp, MatrixBase::operator+, class VectorwiseOp, MatrixBase::sum()
*/
template<typename Scalar> struct ei_scalar_sum_op EIGEN_EMPTY_STRUCT {
EIGEN_STRONG_INLINE const Scalar operator() (const Scalar& a, const Scalar& b) const { return a + b; }
@ -52,7 +52,7 @@ struct ei_functor_traits<ei_scalar_sum_op<Scalar> > {
/** \internal
* \brief Template functor to compute the product of two scalars
*
* \sa class CwiseBinaryOp, Cwise::operator*(), class PartialRedux, MatrixBase::redux()
* \sa class CwiseBinaryOp, Cwise::operator*(), class VectorwiseOp, MatrixBase::redux()
*/
template<typename Scalar> struct ei_scalar_product_op EIGEN_EMPTY_STRUCT {
EIGEN_STRONG_INLINE const Scalar operator() (const Scalar& a, const Scalar& b) const { return a * b; }
@ -74,7 +74,7 @@ struct ei_functor_traits<ei_scalar_product_op<Scalar> > {
/** \internal
* \brief Template functor to compute the min of two scalars
*
* \sa class CwiseBinaryOp, MatrixBase::cwiseMin, class PartialRedux, MatrixBase::minCoeff()
* \sa class CwiseBinaryOp, MatrixBase::cwiseMin, class VectorwiseOp, MatrixBase::minCoeff()
*/
template<typename Scalar> struct ei_scalar_min_op EIGEN_EMPTY_STRUCT {
EIGEN_STRONG_INLINE const Scalar operator() (const Scalar& a, const Scalar& b) const { return std::min(a, b); }
@ -96,7 +96,7 @@ struct ei_functor_traits<ei_scalar_min_op<Scalar> > {
/** \internal
* \brief Template functor to compute the max of two scalars
*
* \sa class CwiseBinaryOp, MatrixBase::cwiseMax, class PartialRedux, MatrixBase::maxCoeff()
* \sa class CwiseBinaryOp, MatrixBase::cwiseMax, class VectorwiseOp, MatrixBase::maxCoeff()
*/
template<typename Scalar> struct ei_scalar_max_op EIGEN_EMPTY_STRUCT {
EIGEN_STRONG_INLINE const Scalar operator() (const Scalar& a, const Scalar& b) const { return std::max(a, b); }

View File

@ -32,7 +32,7 @@
*
* This class is the base that is inherited by all matrix, vector, and expression
* types. Most of the Eigen API is contained in this class. Other important classes for
* the Eigen API are Matrix, Cwise, and PartialRedux.
* the Eigen API are Matrix, Cwise, and VectorwiseOp.
*
* Note that some methods are defined in the \ref Array_Module array module.
*
@ -605,8 +605,10 @@ template<typename Derived> class MatrixBase
bool any(void) const;
int count() const;
const PartialRedux<Derived,Horizontal> rowwise() const;
const PartialRedux<Derived,Vertical> colwise() const;
const VectorwiseOp<Derived,Horizontal> rowwise() const;
VectorwiseOp<Derived,Horizontal> rowwise();
const VectorwiseOp<Derived,Vertical> colwise() const;
VectorwiseOp<Derived,Vertical> colwise();
static const CwiseNullaryOp<ei_scalar_random_op<Scalar>,Derived> Random(int rows, int cols);
static const CwiseNullaryOp<ei_scalar_random_op<Scalar>,Derived> Random(int size);

View File

@ -105,7 +105,7 @@ void ei_cache_friendly_product(
// Array module
template<typename ConditionMatrixType, typename ThenMatrixType, typename ElseMatrixType> class Select;
template<typename MatrixType, typename BinaryOp, int Direction> class PartialReduxExpr;
template<typename ExpressionType, int Direction> class PartialRedux;
template<typename ExpressionType, int Direction> class VectorwiseOp;
template<typename MatrixType,int RowFactor,int ColFactor> class Replicate;
template<typename MatrixType, int Direction = BothDirections> class Reverse;

View File

@ -147,13 +147,13 @@ MatrixBase<Derived>::homogeneous() const
* \nonstableyet
* \returns a matrix expression of homogeneous column (or row) vectors
*
* Example: \include PartialRedux_homogeneous.cpp
* Output: \verbinclude PartialRedux_homogeneous.out
* Example: \include VectorwiseOp_homogeneous.cpp
* Output: \verbinclude VectorwiseOp_homogeneous.out
*
* \sa MatrixBase::homogeneous() */
template<typename ExpressionType, int Direction>
inline const Homogeneous<ExpressionType,Direction>
PartialRedux<ExpressionType,Direction>::homogeneous() const
VectorwiseOp<ExpressionType,Direction>::homogeneous() const
{
return _expression();
}
@ -165,7 +165,7 @@ PartialRedux<ExpressionType,Direction>::homogeneous() const
* Example: \include MatrixBase_hnormalized.cpp
* Output: \verbinclude MatrixBase_hnormalized.out
*
* \sa PartialRedux::hnormalized() */
* \sa VectorwiseOp::hnormalized() */
template<typename Derived>
inline const typename MatrixBase<Derived>::HNormalizedReturnType
MatrixBase<Derived>::hnormalized() const
@ -185,8 +185,8 @@ MatrixBase<Derived>::hnormalized() const
*
* \sa MatrixBase::hnormalized() */
template<typename ExpressionType, int Direction>
inline const typename PartialRedux<ExpressionType,Direction>::HNormalizedReturnType
PartialRedux<ExpressionType,Direction>::hnormalized() const
inline const typename VectorwiseOp<ExpressionType,Direction>::HNormalizedReturnType
VectorwiseOp<ExpressionType,Direction>::hnormalized() const
{
return HNormalized_Block(_expression(),0,0,
Direction==Vertical ? _expression().rows()-1 : _expression().rows(),

View File

@ -105,8 +105,8 @@ MatrixBase<Derived>::cross3(const MatrixBase<OtherDerived>& other) const
* \sa MatrixBase::cross() */
template<typename ExpressionType, int Direction>
template<typename OtherDerived>
const typename PartialRedux<ExpressionType,Direction>::CrossReturnType
PartialRedux<ExpressionType,Direction>::cross(const MatrixBase<OtherDerived>& other) const
const typename VectorwiseOp<ExpressionType,Direction>::CrossReturnType
VectorwiseOp<ExpressionType,Direction>::cross(const MatrixBase<OtherDerived>& other) const
{
EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3)
EIGEN_STATIC_ASSERT((ei_is_same_type<Scalar, typename OtherDerived::Scalar>::ret),

View File

@ -205,7 +205,7 @@ public:
* If Mode==Affine, then the last row is set to [0 ... 0 1] */
inline Transform()
{
if (Mode==Affine)
if (int(Mode)==Affine)
makeAffine();
}

View File

@ -562,8 +562,8 @@ template<typename Derived> class SparseMatrixBase
bool all(void) const;
bool any(void) const;
const PartialRedux<Derived,Horizontal> rowwise() const;
const PartialRedux<Derived,Vertical> colwise() const;
const VectorwiseOp<Derived,Horizontal> rowwise() const;
const VectorwiseOp<Derived,Vertical> colwise() const;
static const CwiseNullaryOp<ei_scalar_random_op<Scalar>,Derived> Random(int rows, int cols);
static const CwiseNullaryOp<ei_scalar_random_op<Scalar>,Derived> Random(int size);

View File

@ -33,7 +33,8 @@ template<typename MatrixType> void array(const MatrixType& m)
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> ColVectorType;
typedef Matrix<Scalar, 1, MatrixType::ColsAtCompileTime> RowVectorType;
int rows = m.rows();
int cols = m.cols();
@ -42,6 +43,9 @@ template<typename MatrixType> void array(const MatrixType& m)
m2 = MatrixType::Random(rows, cols),
m3(rows, cols);
ColVectorType cv1 = ColVectorType::Random(rows);
RowVectorType rv1 = RowVectorType::Random(cols);
Scalar s1 = ei_random<Scalar>(),
s2 = ei_random<Scalar>();
@ -62,6 +66,16 @@ template<typename MatrixType> void array(const MatrixType& m)
if (!ei_isApprox(m1.sum(), (m1+m2).sum()))
VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(ei_scalar_sum_op<Scalar>()));
// vector-wise ops
m3 = m1;
VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
m3 = m1;
VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
m3 = m1;
VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
m3 = m1;
VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
}
template<typename MatrixType> void comparisons(const MatrixType& m)