adapt CwiseUnaryOp and CwiseUnaryView

This commit is contained in:
Gael Guennebaud 2009-11-16 19:39:29 +01:00
parent 2a3a6fe45e
commit 1c9a2d246f
6 changed files with 341 additions and 311 deletions

View File

@ -56,13 +56,17 @@ struct ei_traits<CwiseUnaryOp<UnaryOp, MatrixType> >
};
};
template<typename UnaryOp, typename MatrixType, typename StorageType>
class CwiseUnaryOpImpl;
template<typename UnaryOp, typename MatrixType>
class CwiseUnaryOp : ei_no_assignment_operator,
public MatrixBase<CwiseUnaryOp<UnaryOp, MatrixType> >
public CwiseUnaryOpImpl<UnaryOp, MatrixType, typename ei_traits<MatrixType>::StorageType>
{
public:
EIGEN_GENERIC_PUBLIC_INTERFACE(CwiseUnaryOp)
typedef typename CwiseUnaryOpImpl<UnaryOp, MatrixType,typename ei_traits<MatrixType>::StorageType>::Base Base;
EIGEN_GENERIC_PUBLIC_INTERFACE_NEW(CwiseUnaryOp)
inline CwiseUnaryOp(const MatrixType& mat, const UnaryOp& func = UnaryOp())
: m_matrix(mat), m_functor(func) {}
@ -70,28 +74,6 @@ class CwiseUnaryOp : ei_no_assignment_operator,
EIGEN_STRONG_INLINE int rows() const { return m_matrix.rows(); }
EIGEN_STRONG_INLINE int cols() const { return m_matrix.cols(); }
EIGEN_STRONG_INLINE const Scalar coeff(int row, int col) const
{
return m_functor(m_matrix.coeff(row, col));
}
template<int LoadMode>
EIGEN_STRONG_INLINE PacketScalar packet(int row, int col) const
{
return m_functor.packetOp(m_matrix.template packet<LoadMode>(row, col));
}
EIGEN_STRONG_INLINE const Scalar coeff(int index) const
{
return m_functor(m_matrix.coeff(index));
}
template<int LoadMode>
EIGEN_STRONG_INLINE PacketScalar packet(int index) const
{
return m_functor.packetOp(m_matrix.template packet<LoadMode>(index));
}
/** \internal used for introspection */
const UnaryOp& _functor() const { return m_functor; }
@ -99,39 +81,53 @@ class CwiseUnaryOp : ei_no_assignment_operator,
const typename ei_cleantype<typename MatrixType::Nested>::type&
_expression() const { return m_matrix; }
const typename ei_cleantype<typename MatrixType::Nested>::type&
nestedExpression() const { return m_matrix; }
typename ei_cleantype<typename MatrixType::Nested>::type&
nestedExpression() { return m_matrix.const_cast_derived(); }
protected:
const typename MatrixType::Nested m_matrix;
const UnaryOp m_functor;
};
/** \returns an expression of a custom coefficient-wise unary operator \a func of *this
*
* The template parameter \a CustomUnaryOp is the type of the functor
* of the custom unary operator.
*
* Example:
* \include class_CwiseUnaryOp.cpp
* Output: \verbinclude class_CwiseUnaryOp.out
*
* \sa class CwiseUnaryOp, class CwiseBinarOp, MatrixBase::operator-, Cwise::abs
*/
template<typename Derived>
template<typename CustomUnaryOp>
EIGEN_STRONG_INLINE const CwiseUnaryOp<CustomUnaryOp, Derived>
MatrixBase<Derived>::unaryExpr(const CustomUnaryOp& func) const
template<typename UnaryOp, typename MatrixType>
class CwiseUnaryOpImpl<UnaryOp,MatrixType,Dense> : public MatrixBase<CwiseUnaryOp<UnaryOp, MatrixType> >
{
return CwiseUnaryOp<CustomUnaryOp, Derived>(derived(), func);
const typename ei_cleantype<typename MatrixType::Nested>::type& matrix() const
{ return derived().nestedExpression(); }
typename ei_cleantype<typename MatrixType::Nested>::type& matrix()
{ return derived().nestedExpression(); }
public:
typedef CwiseUnaryOp<UnaryOp, MatrixType> Derived;
EIGEN_DENSE_PUBLIC_INTERFACE( Derived )
EIGEN_STRONG_INLINE const Scalar coeff(int row, int col) const
{
return derived()._functor()(matrix().coeff(row, col));
}
/** \returns an expression of the opposite of \c *this
*/
template<typename Derived>
EIGEN_STRONG_INLINE const CwiseUnaryOp<ei_scalar_opposite_op<typename ei_traits<Derived>::Scalar>,Derived>
MatrixBase<Derived>::operator-() const
template<int LoadMode>
EIGEN_STRONG_INLINE PacketScalar packet(int row, int col) const
{
return derived();
return derived()._functor().packetOp(matrix().template packet<LoadMode>(row, col));
}
EIGEN_STRONG_INLINE const Scalar coeff(int index) const
{
return derived()._functor()(matrix().coeff(index));
}
template<int LoadMode>
EIGEN_STRONG_INLINE PacketScalar packet(int index) const
{
return derived()._functor().packetOp(matrix().template packet<LoadMode>(index));
}
};
/** \returns an expression of the coefficient-wise absolute value of \c *this
*
* Example: \include Cwise_abs.cpp
@ -160,49 +156,6 @@ Cwise<ExpressionType>::abs2() const
return _expression();
}
/** \returns an expression of the complex conjugate of \c *this.
*
* \sa adjoint() */
template<typename Derived>
EIGEN_STRONG_INLINE typename MatrixBase<Derived>::ConjugateReturnType
MatrixBase<Derived>::conjugate() const
{
return ConjugateReturnType(derived());
}
/** \returns a read-only expression of the real part of \c *this.
*
* \sa imag() */
template<typename Derived>
EIGEN_STRONG_INLINE typename MatrixBase<Derived>::RealReturnType
MatrixBase<Derived>::real() const { return derived(); }
/** \returns an read-only expression of the imaginary part of \c *this.
*
* \sa real() */
template<typename Derived>
EIGEN_STRONG_INLINE const typename MatrixBase<Derived>::ImagReturnType
MatrixBase<Derived>::imag() const { return derived(); }
/** \returns an expression of *this with the \a Scalar type casted to
* \a NewScalar.
*
* The template parameter \a NewScalar is the type we are casting the scalars to.
*
* \sa class CwiseUnaryOp
*/
template<typename Derived>
template<typename NewType>
EIGEN_STRONG_INLINE
typename ei_cast_return_type<
Derived,
const CwiseUnaryOp<ei_scalar_cast_op<typename ei_traits<Derived>::Scalar, NewType>, Derived>
>::type
MatrixBase<Derived>::cast() const
{
return derived();
}
/** \returns an expression of the coefficient-wise exponential of *this.
*
* Example: \include Cwise_exp.cpp
@ -231,47 +184,4 @@ Cwise<ExpressionType>::log() const
return _expression();
}
/** \returns an expression of \c *this scaled by the scalar factor \a scalar */
template<typename Derived>
EIGEN_STRONG_INLINE const typename MatrixBase<Derived>::ScalarMultipleReturnType
MatrixBase<Derived>::operator*(const Scalar& scalar) const
{
return CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, Derived>
(derived(), ei_scalar_multiple_op<Scalar>(scalar));
}
/** Overloaded for efficient real matrix times complex scalar value */
template<typename Derived>
EIGEN_STRONG_INLINE const CwiseUnaryOp<ei_scalar_multiple2_op<typename ei_traits<Derived>::Scalar,
std::complex<typename ei_traits<Derived>::Scalar> >, Derived>
MatrixBase<Derived>::operator*(const std::complex<Scalar>& scalar) const
{
return CwiseUnaryOp<ei_scalar_multiple2_op<Scalar,std::complex<Scalar> >, Derived>
(*static_cast<const Derived*>(this), ei_scalar_multiple2_op<Scalar,std::complex<Scalar> >(scalar));
}
/** \returns an expression of \c *this divided by the scalar value \a scalar */
template<typename Derived>
EIGEN_STRONG_INLINE const CwiseUnaryOp<ei_scalar_quotient1_op<typename ei_traits<Derived>::Scalar>, Derived>
MatrixBase<Derived>::operator/(const Scalar& scalar) const
{
return CwiseUnaryOp<ei_scalar_quotient1_op<Scalar>, Derived>
(derived(), ei_scalar_quotient1_op<Scalar>(scalar));
}
template<typename Derived>
EIGEN_STRONG_INLINE Derived&
MatrixBase<Derived>::operator*=(const Scalar& other)
{
return *this = *this * other;
}
template<typename Derived>
EIGEN_STRONG_INLINE Derived&
MatrixBase<Derived>::operator/=(const Scalar& other)
{
return *this = *this / other;
}
#endif // EIGEN_CWISE_UNARY_OP_H

View File

@ -0,0 +1,158 @@
#ifndef EIGEN_PARSED_BY_DOXYGEN
/** \internal Represents a scalar multiple of a matrix */
typedef CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, Derived> ScalarMultipleReturnType;
/** \internal Represents a quotient of a matrix by a scalar*/
typedef CwiseUnaryOp<ei_scalar_quotient1_op<Scalar>, Derived> ScalarQuotient1ReturnType;
/** \internal the return type of MatrixBase::conjugate() */
typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
const CwiseUnaryOp<ei_scalar_conjugate_op<Scalar>, Derived>,
const Derived&
>::ret ConjugateReturnType;
/** \internal the return type of MatrixBase::real() const */
typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
const CwiseUnaryOp<ei_scalar_real_op<Scalar>, Derived>,
const Derived&
>::ret RealReturnType;
/** \internal the return type of MatrixBase::real() */
typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
CwiseUnaryView<ei_scalar_real_op<Scalar>, Derived>,
Derived&
>::ret NonConstRealReturnType;
/** \internal the return type of MatrixBase::imag() const */
typedef CwiseUnaryOp<ei_scalar_imag_op<Scalar>, Derived> ImagReturnType;
/** \internal the return type of MatrixBase::imag() */
typedef CwiseUnaryView<ei_scalar_imag_op<Scalar>, Derived> NonConstImagReturnType;
#endif // not EIGEN_PARSED_BY_DOXYGEN
/** \returns an expression of the opposite of \c *this
*/
EIGEN_STRONG_INLINE const CwiseUnaryOp<ei_scalar_opposite_op<typename ei_traits<Derived>::Scalar>,Derived>
operator-() const { return derived(); }
EIGEN_STRONG_INLINE Derived& operator*=(const Scalar& other)
{ return *this = *this * other; }
EIGEN_STRONG_INLINE Derived& operator/=(const Scalar& other)
{ return *this = *this / other; }
/** \returns an expression of \c *this scaled by the scalar factor \a scalar */
EIGEN_STRONG_INLINE const ScalarMultipleReturnType
operator*(const Scalar& scalar) const
{
return CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, Derived>
(derived(), ei_scalar_multiple_op<Scalar>(scalar));
}
#ifdef EIGEN_PARSED_BY_DOXYGEN
const ScalarMultipleReturnType operator*(const RealScalar& scalar) const;
#endif
/** \returns an expression of \c *this divided by the scalar value \a scalar */
EIGEN_STRONG_INLINE const CwiseUnaryOp<ei_scalar_quotient1_op<typename ei_traits<Derived>::Scalar>, Derived>
operator/(const Scalar& scalar) const
{
return CwiseUnaryOp<ei_scalar_quotient1_op<Scalar>, Derived>
(derived(), ei_scalar_quotient1_op<Scalar>(scalar));
}
/** Overloaded for efficient real matrix times complex scalar value */
EIGEN_STRONG_INLINE const CwiseUnaryOp<ei_scalar_multiple2_op<Scalar,std::complex<Scalar> >, Derived>
operator*(const std::complex<Scalar>& scalar) const
{
return CwiseUnaryOp<ei_scalar_multiple2_op<Scalar,std::complex<Scalar> >, Derived>
(*static_cast<const Derived*>(this), ei_scalar_multiple2_op<Scalar,std::complex<Scalar> >(scalar));
}
inline friend const ScalarMultipleReturnType
operator*(const Scalar& scalar, const MatrixBase& matrix)
{ return matrix*scalar; }
inline friend const CwiseUnaryOp<ei_scalar_multiple2_op<Scalar,std::complex<Scalar> >, Derived>
operator*(const std::complex<Scalar>& scalar, const MatrixBase& matrix)
{ return matrix*scalar; }
/** \returns an expression of *this with the \a Scalar type casted to
* \a NewScalar.
*
* The template parameter \a NewScalar is the type we are casting the scalars to.
*
* \sa class CwiseUnaryOp
*/
template<typename NewType>
typename ei_cast_return_type<Derived,const CwiseUnaryOp<ei_scalar_cast_op<typename ei_traits<Derived>::Scalar, NewType>, Derived> >::type
cast() const
{
return derived();
}
/** \returns an expression of the complex conjugate of \c *this.
*
* \sa adjoint() */
EIGEN_STRONG_INLINE ConjugateReturnType
conjugate() const
{
return ConjugateReturnType(derived());
}
/** \returns a read-only expression of the real part of \c *this.
*
* \sa imag() */
EIGEN_STRONG_INLINE RealReturnType
real() const { return derived(); }
/** \returns an read-only expression of the imaginary part of \c *this.
*
* \sa real() */
EIGEN_STRONG_INLINE const ImagReturnType
imag() const { return derived(); }
/** \returns an expression of a custom coefficient-wise unary operator \a func of *this
*
* The template parameter \a CustomUnaryOp is the type of the functor
* of the custom unary operator.
*
* Example:
* \include class_CwiseUnaryOp.cpp
* Output: \verbinclude class_CwiseUnaryOp.out
*
* \sa class CwiseUnaryOp, class CwiseBinarOp, MatrixBase::operator-, Cwise::abs
*/
template<typename CustomUnaryOp>
EIGEN_STRONG_INLINE const CwiseUnaryOp<CustomUnaryOp, Derived>
unaryExpr(const CustomUnaryOp& func = CustomUnaryOp()) const
{
return CwiseUnaryOp<CustomUnaryOp, Derived>(derived(), func);
}
/** \returns an expression of a custom coefficient-wise unary operator \a func of *this
*
* The template parameter \a CustomUnaryOp is the type of the functor
* of the custom unary operator.
*
* Example:
* \include class_CwiseUnaryOp.cpp
* Output: \verbinclude class_CwiseUnaryOp.out
*
* \sa class CwiseUnaryOp, class CwiseBinarOp, MatrixBase::operator-, Cwise::abs
*/
template<typename CustomViewOp>
EIGEN_STRONG_INLINE const CwiseUnaryView<CustomViewOp, Derived>
unaryViewExpr(const CustomViewOp& func = CustomViewOp()) const
{
return CwiseUnaryView<CustomViewOp, Derived>(derived(), func);
}
/** \returns a non const expression of the real part of \c *this.
*
* \sa imag() */
EIGEN_STRONG_INLINE NonConstRealReturnType
real() { return derived(); }
/** \returns a non const expression of the imaginary part of \c *this.
*
* \sa real() */
EIGEN_STRONG_INLINE NonConstImagReturnType
imag() { return derived(); }
const Cwise<Derived> cwise() const;
Cwise<Derived> cwise();

View File

@ -52,12 +52,17 @@ struct ei_traits<CwiseUnaryView<ViewOp, MatrixType> >
};
};
template<typename ViewOp, typename MatrixType, typename StorageType>
class CwiseUnaryViewImpl;
template<typename ViewOp, typename MatrixType>
class CwiseUnaryView : public MatrixBase<CwiseUnaryView<ViewOp, MatrixType> >
class CwiseUnaryView : ei_no_assignment_operator,
public CwiseUnaryViewImpl<ViewOp, MatrixType, typename ei_traits<MatrixType>::StorageType>
{
public:
EIGEN_GENERIC_PUBLIC_INTERFACE(CwiseUnaryView)
typedef typename CwiseUnaryViewImpl<ViewOp, MatrixType,typename ei_traits<MatrixType>::StorageType>::Base Base;
EIGEN_GENERIC_PUBLIC_INTERFACE_NEW(CwiseUnaryView)
inline CwiseUnaryView(const MatrixType& mat, const ViewOp& func = ViewOp())
: m_matrix(mat), m_functor(func) {}
@ -67,25 +72,14 @@ class CwiseUnaryView : public MatrixBase<CwiseUnaryView<ViewOp, MatrixType> >
EIGEN_STRONG_INLINE int rows() const { return m_matrix.rows(); }
EIGEN_STRONG_INLINE int cols() const { return m_matrix.cols(); }
EIGEN_STRONG_INLINE const Scalar coeff(int row, int col) const
{
return m_functor(m_matrix.coeff(row, col));
}
/** \internal used for introspection */
const ViewOp& _functor() const { return m_functor; }
EIGEN_STRONG_INLINE const Scalar coeff(int index) const
{
return m_functor(m_matrix.coeff(index));
}
const typename ei_cleantype<typename MatrixType::Nested>::type&
nestedExpression() const { return m_matrix; }
EIGEN_STRONG_INLINE Scalar& coeffRef(int row, int col)
{
return m_functor(m_matrix.const_cast_derived().coeffRef(row, col));
}
EIGEN_STRONG_INLINE Scalar& coeffRef(int index)
{
return m_functor(m_matrix.const_cast_derived().coeffRef(index));
}
typename ei_cleantype<typename MatrixType::Nested>::type&
nestedExpression() { return m_matrix.const_cast_derived(); }
protected:
// FIXME changed from MatrixType::Nested because of a weird compilation error with sun CC
@ -93,37 +87,40 @@ class CwiseUnaryView : public MatrixBase<CwiseUnaryView<ViewOp, MatrixType> >
const ViewOp m_functor;
};
/** \returns an expression of a custom coefficient-wise unary operator \a func of *this
*
* The template parameter \a CustomUnaryOp is the type of the functor
* of the custom unary operator.
*
* Example:
* \include class_CwiseUnaryOp.cpp
* Output: \verbinclude class_CwiseUnaryOp.out
*
* \sa class CwiseUnaryOp, class CwiseBinarOp, MatrixBase::operator-, Cwise::abs
*/
template<typename Derived>
template<typename CustomViewOp>
EIGEN_STRONG_INLINE const CwiseUnaryView<CustomViewOp, Derived>
MatrixBase<Derived>::unaryViewExpr(const CustomViewOp& func) const
template<typename ViewOp, typename MatrixType>
class CwiseUnaryViewImpl<ViewOp,MatrixType,Dense> : public MatrixBase<CwiseUnaryView<ViewOp, MatrixType> >
{
return CwiseUnaryView<CustomViewOp, Derived>(derived(), func);
const typename ei_cleantype<typename MatrixType::Nested>::type& matrix() const
{ return derived().nestedExpression(); }
typename ei_cleantype<typename MatrixType::Nested>::type& matrix()
{ return derived().nestedExpression(); }
public:
typedef CwiseUnaryView<ViewOp, MatrixType> Derived;
EIGEN_DENSE_PUBLIC_INTERFACE( Derived )
EIGEN_STRONG_INLINE const Scalar coeff(int row, int col) const
{
return derived()._functor()(matrix().coeff(row, col));
}
/** \returns a non const expression of the real part of \c *this.
*
* \sa imag() */
template<typename Derived>
EIGEN_STRONG_INLINE typename MatrixBase<Derived>::NonConstRealReturnType
MatrixBase<Derived>::real() { return derived(); }
EIGEN_STRONG_INLINE const Scalar coeff(int index) const
{
return derived()._functor()(matrix().coeff(index));
}
EIGEN_STRONG_INLINE Scalar& coeffRef(int row, int col)
{
return derived()._functor()(matrix().const_cast_derived().coeffRef(row, col));
}
EIGEN_STRONG_INLINE Scalar& coeffRef(int index)
{
return derived()._functor()(matrix().const_cast_derived().coeffRef(index));
}
};
/** \returns a non const expression of the imaginary part of \c *this.
*
* \sa real() */
template<typename Derived>
EIGEN_STRONG_INLINE typename MatrixBase<Derived>::NonConstImagReturnType
MatrixBase<Derived>::imag() { return derived(); }
#endif // EIGEN_CWISE_UNARY_VIEW_H

View File

@ -226,29 +226,6 @@ template<typename Derived> class MatrixBase
/** \internal Represents a matrix with all coefficients equal to one another*/
typedef CwiseNullaryOp<ei_scalar_constant_op<Scalar>,Derived> ConstantReturnType;
/** \internal Represents a scalar multiple of a matrix */
typedef CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, Derived> ScalarMultipleReturnType;
/** \internal Represents a quotient of a matrix by a scalar*/
typedef CwiseUnaryOp<ei_scalar_quotient1_op<Scalar>, Derived> ScalarQuotient1ReturnType;
/** \internal the return type of MatrixBase::conjugate() */
typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
const CwiseUnaryOp<ei_scalar_conjugate_op<Scalar>, Derived>,
const Derived&
>::ret ConjugateReturnType;
/** \internal the return type of MatrixBase::real() const */
typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
const CwiseUnaryOp<ei_scalar_real_op<Scalar>, Derived>,
const Derived&
>::ret RealReturnType;
/** \internal the return type of MatrixBase::real() */
typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
CwiseUnaryView<ei_scalar_real_op<Scalar>, Derived>,
Derived&
>::ret NonConstRealReturnType;
/** \internal the return type of MatrixBase::imag() const */
typedef CwiseUnaryOp<ei_scalar_imag_op<Scalar>, Derived> ImagReturnType;
/** \internal the return type of MatrixBase::imag() */
typedef CwiseUnaryView<ei_scalar_imag_op<Scalar>, Derived> NonConstImagReturnType;
/** \internal the return type of MatrixBase::adjoint() */
typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
CwiseUnaryOp<ei_scalar_conjugate_op<Scalar>, NestByValue<Eigen::Transpose<Derived> > >,
@ -268,6 +245,7 @@ template<typename Derived> class MatrixBase
ei_traits<Derived>::ColsAtCompileTime> BasisReturnType;
#endif // not EIGEN_PARSED_BY_DOXYGEN
#include "CwiseUnaryOps.h"
/** Copies \a other into *this. \returns a reference to *this. */
template<typename OtherDerived>
@ -363,8 +341,6 @@ template<typename Derived> class MatrixBase
Scalar& w();
const CwiseUnaryOp<ei_scalar_opposite_op<typename ei_traits<Derived>::Scalar>,Derived> operator-() const;
template<typename OtherDerived>
const CwiseBinaryOp<ei_scalar_sum_op<typename ei_traits<Derived>::Scalar>, Derived, OtherDerived>
operator+(const MatrixBase<OtherDerived> &other) const;
@ -378,27 +354,6 @@ template<typename Derived> class MatrixBase
template<typename OtherDerived>
Derived& operator-=(const MatrixBase<OtherDerived>& other);
Derived& operator*=(const Scalar& other);
Derived& operator/=(const Scalar& other);
const ScalarMultipleReturnType operator*(const Scalar& scalar) const;
#ifdef EIGEN_PARSED_BY_DOXYGEN
const ScalarMultipleReturnType operator*(const RealScalar& scalar) const;
#endif
const CwiseUnaryOp<ei_scalar_quotient1_op<typename ei_traits<Derived>::Scalar>, Derived>
operator/(const Scalar& scalar) const;
const CwiseUnaryOp<ei_scalar_multiple2_op<Scalar,std::complex<Scalar> >, Derived>
operator*(const std::complex<Scalar>& scalar) const;
inline friend const ScalarMultipleReturnType
operator*(const Scalar& scalar, const MatrixBase& matrix)
{ return matrix*scalar; }
inline friend const CwiseUnaryOp<ei_scalar_multiple2_op<Scalar,std::complex<Scalar> >, Derived>
operator*(const std::complex<Scalar>& scalar, const MatrixBase& matrix)
{ return matrix*scalar; }
template<typename OtherDerived>
const typename ProductReturnType<Derived,OtherDerived>::Type
operator*(const MatrixBase<OtherDerived> &other) const;
@ -592,13 +547,6 @@ template<typename Derived> class MatrixBase
{ return (cwise() != other).any(); }
template<typename NewType>
typename ei_cast_return_type<
Derived,
const CwiseUnaryOp<ei_scalar_cast_op<typename ei_traits<Derived>::Scalar, NewType>, Derived>
>::type
cast() const;
/** \returns the matrix or vector obtained by evaluating this expression.
*
* Notice that in the case of a plain matrix or vector (not an expression) this function just returns
@ -626,18 +574,6 @@ template<typename Derived> class MatrixBase
inline const NestByValue<Derived> nestByValue() const;
ConjugateReturnType conjugate() const;
RealReturnType real() const;
NonConstRealReturnType real();
const ImagReturnType imag() const;
NonConstImagReturnType imag();
template<typename CustomUnaryOp>
const CwiseUnaryOp<CustomUnaryOp, Derived> unaryExpr(const CustomUnaryOp& func = CustomUnaryOp()) const;
template<typename CustomViewOp>
const CwiseUnaryView<CustomViewOp, Derived> unaryViewExpr(const CustomViewOp& func = CustomViewOp()) const;
template<typename CustomBinaryOp, typename OtherDerived>
const CwiseBinaryOp<CustomBinaryOp, Derived, OtherDerived>
binaryExpr(const MatrixBase<OtherDerived> &other, const CustomBinaryOp& func = CustomBinaryOp()) const;
@ -671,9 +607,6 @@ template<typename Derived> class MatrixBase
{ return *static_cast<Derived*>(const_cast<MatrixBase*>(this)); }
#endif // not EIGEN_PARSED_BY_DOXYGEN
const Cwise<Derived> cwise() const;
Cwise<Derived> cwise();
inline const WithFormat<Derived> format(const IOFormat& fmt) const;
/////////// Array module ///////////

View File

@ -25,55 +25,42 @@
#ifndef EIGEN_SPARSE_CWISE_UNARY_OP_H
#define EIGEN_SPARSE_CWISE_UNARY_OP_H
template<typename UnaryOp, typename MatrixType>
struct ei_traits<SparseCwiseUnaryOp<UnaryOp, MatrixType> > : ei_traits<MatrixType>
{
typedef typename ei_result_of<
UnaryOp(typename MatrixType::Scalar)
>::type Scalar;
typedef typename MatrixType::Nested MatrixTypeNested;
typedef typename ei_unref<MatrixTypeNested>::type _MatrixTypeNested;
enum {
CoeffReadCost = _MatrixTypeNested::CoeffReadCost + ei_functor_traits<UnaryOp>::Cost
};
};
// template<typename UnaryOp, typename MatrixType>
// struct ei_traits<SparseCwiseUnaryOp<UnaryOp, MatrixType> > : ei_traits<MatrixType>
// {
// typedef typename ei_result_of<
// UnaryOp(typename MatrixType::Scalar)
// >::type Scalar;
// typedef typename MatrixType::Nested MatrixTypeNested;
// typedef typename ei_unref<MatrixTypeNested>::type _MatrixTypeNested;
// enum {
// CoeffReadCost = _MatrixTypeNested::CoeffReadCost + ei_functor_traits<UnaryOp>::Cost
// };
// };
template<typename UnaryOp, typename MatrixType>
class SparseCwiseUnaryOp : ei_no_assignment_operator,
public SparseMatrixBase<SparseCwiseUnaryOp<UnaryOp, MatrixType> >
class CwiseUnaryOpImp<UnaryOp,MatrixType,Sparse>
: public SparseMatrixBase<CwiseUnaryOp<UnaryOp, MatrixType> >
{
public:
class InnerIterator;
// typedef typename ei_unref<LhsNested>::type _LhsNested;
EIGEN_SPARSE_GENERIC_PUBLIC_INTERFACE(SparseCwiseUnaryOp)
inline SparseCwiseUnaryOp(const MatrixType& mat, const UnaryOp& func = UnaryOp())
: m_matrix(mat), m_functor(func) {}
EIGEN_STRONG_INLINE int rows() const { return m_matrix.rows(); }
EIGEN_STRONG_INLINE int cols() const { return m_matrix.cols(); }
// EIGEN_STRONG_INLINE const typename MatrixType::Nested& _matrix() const { return m_matrix; }
// EIGEN_STRONG_INLINE const UnaryOp& _functor() const { return m_functor; }
protected:
const typename MatrixType::Nested m_matrix;
const UnaryOp m_functor;
typedef CwiseUnaryOp<UnaryOp, MatrixType> Derived;
EIGEN_SPARSE_PUBLIC_INTERFACE(Derived)
};
template<typename UnaryOp, typename MatrixType>
class SparseCwiseUnaryOp<UnaryOp,MatrixType>::InnerIterator
class CwiseUnaryOpImpl<UnaryOp,MatrixType,Sparse>::InnerIterator
{
typedef typename SparseCwiseUnaryOp::Scalar Scalar;
typedef typename ei_traits<SparseCwiseUnaryOp>::_MatrixTypeNested _MatrixTypeNested;
typedef typename CwiseUnaryOpImpl::Scalar Scalar;
typedef typename ei_traits<CwiseUnaryOpImpl>::_MatrixTypeNested _MatrixTypeNested;
typedef typename _MatrixTypeNested::InnerIterator MatrixTypeIterator;
public:
EIGEN_STRONG_INLINE InnerIterator(const SparseCwiseUnaryOp& unaryOp, int outer)
: m_iter(unaryOp.m_matrix,outer), m_functor(unaryOp.m_functor)
: m_iter(unaryOp.nestedExpression(),outer), m_functor(unaryOp._functor())
{}
EIGEN_STRONG_INLINE InnerIterator& operator++()
@ -92,6 +79,49 @@ class SparseCwiseUnaryOp<UnaryOp,MatrixType>::InnerIterator
const UnaryOp m_functor;
};
template<typename ViewOp, typename MatrixType>
class CwiseUnaryOpImp<ViewOp,MatrixType,Sparse>
: public SparseMatrixBase<CwiseUnaryView<ViewOp, MatrixType> >
{
public:
class InnerIterator;
// typedef typename ei_unref<LhsNested>::type _LhsNested;
typedef CwiseUnaryView<ViewOp, MatrixType> Derived;
EIGEN_SPARSE_PUBLIC_INTERFACE(Derived)
};
template<typename ViewOp, typename MatrixType>
class CwiseUnaryViewImpl<ViewOp,MatrixType,Sparse>::InnerIterator
{
typedef typename CwiseUnaryViewImpl::Scalar Scalar;
typedef typename ei_traits<SparseCwiseUnaryOp>::_MatrixTypeNested _MatrixTypeNested;
typedef typename _MatrixTypeNested::InnerIterator MatrixTypeIterator;
public:
EIGEN_STRONG_INLINE InnerIterator(const CwiseUnaryViewImpl& unaryView, int outer)
: m_iter(unaryView.nestedExpression(),outer), m_functor(unaryView._functor())
{}
EIGEN_STRONG_INLINE InnerIterator& operator++()
{ ++m_iter; return *this; }
EIGEN_STRONG_INLINE Scalar value() const { return m_functor(m_iter.value()); }
EIGEN_STRONG_INLINE int index() const { return m_iter.index(); }
EIGEN_STRONG_INLINE int row() const { return m_iter.row(); }
EIGEN_STRONG_INLINE int col() const { return m_iter.col(); }
EIGEN_STRONG_INLINE operator bool() const { return m_iter; }
protected:
MatrixTypeIterator m_iter;
const UnaryOp m_functor;
};
/*
template<typename Derived>
template<typename CustomUnaryOp>
EIGEN_STRONG_INLINE const SparseCwiseUnaryOp<CustomUnaryOp, Derived>
@ -178,6 +208,6 @@ SparseMatrixBase<Derived>::operator/=(const Scalar& other)
for (typename Derived::InnerIterator i(derived(),j); i; ++i)
i.valueRef() /= other;
return derived();
}
}*/
#endif // EIGEN_SPARSE_CWISE_UNARY_OP_H

View File

@ -97,21 +97,23 @@ template<typename Derived> class SparseMatrixBase : public AnyMatrixBase<Derived
#endif
};
/** \internal the return type of MatrixBase::conjugate() */
typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
const SparseCwiseUnaryOp<ei_scalar_conjugate_op<Scalar>, Derived>,
const Derived&
>::ret ConjugateReturnType;
/** \internal the return type of MatrixBase::real() */
typedef SparseCwiseUnaryOp<ei_scalar_real_op<Scalar>, Derived> RealReturnType;
/** \internal the return type of MatrixBase::imag() */
typedef SparseCwiseUnaryOp<ei_scalar_imag_op<Scalar>, Derived> ImagReturnType;
/* \internal the return type of MatrixBase::conjugate() */
// typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
// const SparseCwiseUnaryOp<ei_scalar_conjugate_op<Scalar>, Derived>,
// const Derived&
// >::ret ConjugateReturnType;
/* \internal the return type of MatrixBase::real() */
// typedef SparseCwiseUnaryOp<ei_scalar_real_op<Scalar>, Derived> RealReturnType;
/* \internal the return type of MatrixBase::imag() */
// typedef SparseCwiseUnaryOp<ei_scalar_imag_op<Scalar>, Derived> ImagReturnType;
/** \internal the return type of MatrixBase::adjoint() */
typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
SparseCwiseUnaryOp<ei_scalar_conjugate_op<Scalar>, SparseNestByValue<Eigen::Transpose<Derived> > >,
Transpose<Derived>
>::ret AdjointReturnType;
#include "../Core/CwiseUnaryOps.h"
#ifndef EIGEN_PARSED_BY_DOXYGEN
/** This is the "real scalar" type; if the \a Scalar type is already real numbers
* (e.g. int, float or double) then \a RealScalar is just the same as \a Scalar. If
@ -284,7 +286,7 @@ template<typename Derived> class SparseMatrixBase : public AnyMatrixBase<Derived
return s;
}
const SparseCwiseUnaryOp<ei_scalar_opposite_op<typename ei_traits<Derived>::Scalar>,Derived> operator-() const;
// const SparseCwiseUnaryOp<ei_scalar_opposite_op<typename ei_traits<Derived>::Scalar>,Derived> operator-() const;
template<typename OtherDerived>
const SparseCwiseBinaryOp<ei_scalar_sum_op<typename ei_traits<Derived>::Scalar>, Derived, OtherDerived>
@ -302,17 +304,17 @@ template<typename Derived> class SparseMatrixBase : public AnyMatrixBase<Derived
// template<typename Lhs,typename Rhs>
// Derived& operator+=(const Flagged<Product<Lhs,Rhs,CacheFriendlyProduct>, 0, EvalBeforeNestingBit | EvalBeforeAssigningBit>& other);
Derived& operator*=(const Scalar& other);
Derived& operator/=(const Scalar& other);
// Derived& operator*=(const Scalar& other);
// Derived& operator/=(const Scalar& other);
const SparseCwiseUnaryOp<ei_scalar_multiple_op<typename ei_traits<Derived>::Scalar>, Derived>
operator*(const Scalar& scalar) const;
const SparseCwiseUnaryOp<ei_scalar_quotient1_op<typename ei_traits<Derived>::Scalar>, Derived>
operator/(const Scalar& scalar) const;
// const SparseCwiseUnaryOp<ei_scalar_multiple_op<typename ei_traits<Derived>::Scalar>, Derived>
// operator*(const Scalar& scalar) const;
// const SparseCwiseUnaryOp<ei_scalar_quotient1_op<typename ei_traits<Derived>::Scalar>, Derived>
// operator/(const Scalar& scalar) const;
inline friend const SparseCwiseUnaryOp<ei_scalar_multiple_op<typename ei_traits<Derived>::Scalar>, Derived>
operator*(const Scalar& scalar, const SparseMatrixBase& matrix)
{ return matrix*scalar; }
// inline friend const SparseCwiseUnaryOp<ei_scalar_multiple_op<typename ei_traits<Derived>::Scalar>, Derived>
// operator*(const Scalar& scalar, const SparseMatrixBase& matrix)
// { return matrix*scalar; }
// sparse * sparse
@ -543,12 +545,12 @@ template<typename Derived> class SparseMatrixBase : public AnyMatrixBase<Derived
inline const SparseNestByValue<Derived> nestByValue() const;
ConjugateReturnType conjugate() const;
const RealReturnType real() const;
const ImagReturnType imag() const;
// ConjugateReturnType conjugate() const;
// const RealReturnType real() const;
// const ImagReturnType imag() const;
template<typename CustomUnaryOp>
const SparseCwiseUnaryOp<CustomUnaryOp, Derived> unaryExpr(const CustomUnaryOp& func = CustomUnaryOp()) const;
// template<typename CustomUnaryOp>
// const SparseCwiseUnaryOp<CustomUnaryOp, Derived> unaryExpr(const CustomUnaryOp& func = CustomUnaryOp()) const;
// template<typename CustomBinaryOp, typename OtherDerived>
// const CwiseBinaryOp<CustomBinaryOp, Derived, OtherDerived>
@ -572,8 +574,8 @@ template<typename Derived> class SparseMatrixBase : public AnyMatrixBase<Derived
// void visit(Visitor& func) const;
const SparseCwise<Derived> cwise() const;
SparseCwise<Derived> cwise();
// const SparseCwise<Derived> cwise() const;
// SparseCwise<Derived> cwise();
// inline const WithFormat<Derived> format(const IOFormat& fmt) const;