aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Core
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2009-11-18 18:15:19 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2009-11-18 18:15:19 +0100
commite3d890bc5a89798eff50ff6650292b4fa934f72e (patch)
tree3c7d40332a019dce2773fc6a096d046eeab9fb7a /Eigen/src/Core
parent0529ecfe1b43d40e40755a2d856188d3ded2c14e (diff)
Another big refactoring change:
* add a new Eigen2Support module including Cwise, Flagged, and some other deprecated stuff * add a few cwiseXxx functions * adapt a few modules to use cwiseXxx instead of the .cwise() prefix
Diffstat (limited to 'Eigen/src/Core')
-rw-r--r--Eigen/src/Core/Cwise.h186
-rw-r--r--Eigen/src/Core/CwiseBinaryOp.h139
-rw-r--r--Eigen/src/Core/CwiseBinaryOps.h139
-rw-r--r--Eigen/src/Core/CwiseUnaryOp.h56
-rw-r--r--Eigen/src/Core/CwiseUnaryOps.h60
-rw-r--r--Eigen/src/Core/DiagonalMatrix.h6
-rw-r--r--Eigen/src/Core/Dot.h2
-rw-r--r--Eigen/src/Core/Flagged.h149
-rw-r--r--Eigen/src/Core/Fuzzy.h6
-rw-r--r--Eigen/src/Core/MatrixBase.h65
-rw-r--r--Eigen/src/Core/Product.h6
-rw-r--r--Eigen/src/Core/SolveTriangular.h6
-rw-r--r--Eigen/src/Core/StableNorm.h4
-rw-r--r--Eigen/src/Core/util/ForwardDeclarations.h5
14 files changed, 236 insertions, 593 deletions
diff --git a/Eigen/src/Core/Cwise.h b/Eigen/src/Core/Cwise.h
deleted file mode 100644
index 6cf548bc7..000000000
--- a/Eigen/src/Core/Cwise.h
+++ /dev/null
@@ -1,186 +0,0 @@
-// This file is part of Eigen, a lightweight C++ template library
-// for linear algebra.
-//
-// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
-// Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
-//
-// Eigen is free software; you can redistribute it and/or
-// modify it under the terms of the GNU Lesser General Public
-// License as published by the Free Software Foundation; either
-// version 3 of the License, or (at your option) any later version.
-//
-// Alternatively, you can redistribute it and/or
-// modify it under the terms of the GNU General Public License as
-// published by the Free Software Foundation; either version 2 of
-// the License, or (at your option) any later version.
-//
-// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
-// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
-// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
-// GNU General Public License for more details.
-//
-// You should have received a copy of the GNU Lesser General Public
-// License and a copy of the GNU General Public License along with
-// Eigen. If not, see <http://www.gnu.org/licenses/>.
-
-#ifndef EIGEN_CWISE_H
-#define EIGEN_CWISE_H
-
-/** \internal
- * convenient macro to defined the return type of a cwise binary operation */
-#define EIGEN_CWISE_BINOP_RETURN_TYPE(OP) \
- CwiseBinaryOp<OP<typename ei_traits<ExpressionType>::Scalar>, ExpressionType, OtherDerived>
-
-#define EIGEN_CWISE_PRODUCT_RETURN_TYPE \
- CwiseBinaryOp< \
- ei_scalar_product_op< \
- typename ei_scalar_product_traits< \
- typename ei_traits<ExpressionType>::Scalar, \
- typename ei_traits<OtherDerived>::Scalar \
- >::ReturnType \
- >, \
- ExpressionType, \
- OtherDerived \
- >
-
-/** \internal
- * convenient macro to defined the return type of a cwise unary operation */
-#define EIGEN_CWISE_UNOP_RETURN_TYPE(OP) \
- CwiseUnaryOp<OP<typename ei_traits<ExpressionType>::Scalar>, ExpressionType>
-
-/** \internal
- * convenient macro to defined the return type of a cwise comparison to a scalar */
-#define EIGEN_CWISE_COMP_TO_SCALAR_RETURN_TYPE(OP) \
- CwiseBinaryOp<OP<typename ei_traits<ExpressionType>::Scalar>, ExpressionType, \
- NestByValue<typename ExpressionType::ConstantReturnType> >
-
-/** \class Cwise
- *
- * \brief Pseudo expression providing additional coefficient-wise operations
- *
- * \param ExpressionType the type of the object on which to do coefficient-wise operations
- *
- * This class represents an expression with additional coefficient-wise features.
- * It is the return type of MatrixBase::cwise()
- * and most of the time this is the only way it is used.
- *
- * Note that some methods are defined in the \ref Array_Module array module.
- *
- * Example: \include MatrixBase_cwise_const.cpp
- * Output: \verbinclude MatrixBase_cwise_const.out
- *
- * \sa MatrixBase::cwise() const, MatrixBase::cwise()
- */
-template<typename ExpressionType, template<typename> class StorageBase> class Cwise
-{
- public:
-
- typedef typename ei_traits<ExpressionType>::Scalar Scalar;
- typedef typename ei_meta_if<ei_must_nest_by_value<ExpressionType>::ret,
- ExpressionType, const ExpressionType&>::ret ExpressionTypeNested;
- typedef CwiseUnaryOp<ei_scalar_add_op<Scalar>, ExpressionType> ScalarAddReturnType;
-
- inline Cwise(const ExpressionType& matrix) : m_matrix(matrix) {}
-
- /** \internal */
- inline const ExpressionType& _expression() const { return m_matrix; }
-
- template<typename OtherDerived>
- const EIGEN_CWISE_PRODUCT_RETURN_TYPE
- operator*(const AnyMatrixBase<OtherDerived> &other) const;
-
- template<typename OtherDerived>
- const EIGEN_CWISE_BINOP_RETURN_TYPE(ei_scalar_quotient_op)
- operator/(const StorageBase<OtherDerived> &other) const;
-
- template<typename OtherDerived>
- const EIGEN_CWISE_BINOP_RETURN_TYPE(ei_scalar_min_op)
- min(const StorageBase<OtherDerived> &other) const;
-
- template<typename OtherDerived>
- const EIGEN_CWISE_BINOP_RETURN_TYPE(ei_scalar_max_op)
- max(const StorageBase<OtherDerived> &other) const;
-
- const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_abs_op) abs() const;
- const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_abs2_op) abs2() const;
- const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_square_op) square() const;
- const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_cube_op) cube() const;
- const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_inverse_op) inverse() const;
- const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_sqrt_op) sqrt() const;
- const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_exp_op) exp() const;
- const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_log_op) log() const;
- const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_cos_op) cos() const;
- const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_sin_op) sin() const;
- const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_pow_op) pow(const Scalar& exponent) const;
-
- const ScalarAddReturnType
- operator+(const Scalar& scalar) const;
-
- /** \relates Cwise */
- friend const ScalarAddReturnType
- operator+(const Scalar& scalar, const Cwise& mat)
- { return mat + scalar; }
-
- ExpressionType& operator+=(const Scalar& scalar);
-
- const ScalarAddReturnType
- operator-(const Scalar& scalar) const;
-
- ExpressionType& operator-=(const Scalar& scalar);
-
- template<typename OtherDerived>
- inline ExpressionType& operator*=(const StorageBase<OtherDerived> &other);
-
- template<typename OtherDerived>
- inline ExpressionType& operator/=(const StorageBase<OtherDerived> &other);
-
- template<typename OtherDerived> const EIGEN_CWISE_BINOP_RETURN_TYPE(std::less)
- operator<(const StorageBase<OtherDerived>& other) const;
-
- template<typename OtherDerived> const EIGEN_CWISE_BINOP_RETURN_TYPE(std::less_equal)
- operator<=(const StorageBase<OtherDerived>& other) const;
-
- template<typename OtherDerived> const EIGEN_CWISE_BINOP_RETURN_TYPE(std::greater)
- operator>(const StorageBase<OtherDerived>& other) const;
-
- template<typename OtherDerived> const EIGEN_CWISE_BINOP_RETURN_TYPE(std::greater_equal)
- operator>=(const StorageBase<OtherDerived>& other) const;
-
- template<typename OtherDerived> const EIGEN_CWISE_BINOP_RETURN_TYPE(std::equal_to)
- operator==(const StorageBase<OtherDerived>& other) const;
-
- template<typename OtherDerived> const EIGEN_CWISE_BINOP_RETURN_TYPE(std::not_equal_to)
- operator!=(const StorageBase<OtherDerived>& other) const;
-
- // comparisons to a scalar value
- const EIGEN_CWISE_COMP_TO_SCALAR_RETURN_TYPE(std::less)
- operator<(Scalar s) const;
-
- const EIGEN_CWISE_COMP_TO_SCALAR_RETURN_TYPE(std::less_equal)
- operator<=(Scalar s) const;
-
- const EIGEN_CWISE_COMP_TO_SCALAR_RETURN_TYPE(std::greater)
- operator>(Scalar s) const;
-
- const EIGEN_CWISE_COMP_TO_SCALAR_RETURN_TYPE(std::greater_equal)
- operator>=(Scalar s) const;
-
- const EIGEN_CWISE_COMP_TO_SCALAR_RETURN_TYPE(std::equal_to)
- operator==(Scalar s) const;
-
- const EIGEN_CWISE_COMP_TO_SCALAR_RETURN_TYPE(std::not_equal_to)
- operator!=(Scalar s) const;
-
- // allow to extend Cwise outside Eigen
- #ifdef EIGEN_CWISE_PLUGIN
- #include EIGEN_CWISE_PLUGIN
- #endif
-
- protected:
- ExpressionTypeNested m_matrix;
-
- private:
- Cwise& operator=(const Cwise&);
-};
-
-#endif // EIGEN_CWISE_H
diff --git a/Eigen/src/Core/CwiseBinaryOp.h b/Eigen/src/Core/CwiseBinaryOp.h
index bba66f2f3..462e0f92d 100644
--- a/Eigen/src/Core/CwiseBinaryOp.h
+++ b/Eigen/src/Core/CwiseBinaryOp.h
@@ -169,22 +169,6 @@ class CwiseBinaryOpImpl<BinaryOp, Lhs, Rhs, Dense>
}
};
-/**\returns an expression of the difference of \c *this and \a other
- *
- * \note If you want to substract a given scalar from all coefficients, see Cwise::operator-().
- *
- * \sa class CwiseBinaryOp, MatrixBase::operator-=(), Cwise::operator-()
- */
-template<typename Derived>
-template<typename OtherDerived>
-EIGEN_STRONG_INLINE const CwiseBinaryOp<ei_scalar_difference_op<typename ei_traits<Derived>::Scalar>,
- Derived, OtherDerived>
-MatrixBase<Derived>::operator-(const MatrixBase<OtherDerived> &other) const
-{
- return CwiseBinaryOp<ei_scalar_difference_op<Scalar>,
- Derived, OtherDerived>(derived(), other.derived());
-}
-
/** replaces \c *this by \c *this - \a other.
*
* \returns a reference to \c *this
@@ -197,22 +181,6 @@ MatrixBase<Derived>::operator-=(const MatrixBase<OtherDerived> &other)
return *this = *this - other;
}
-/** \relates MatrixBase
- *
- * \returns an expression of the sum of \c *this and \a other
- *
- * \note If you want to add a given scalar to all coefficients, see Cwise::operator+().
- *
- * \sa class CwiseBinaryOp, MatrixBase::operator+=(), Cwise::operator+()
- */
-template<typename Derived>
-template<typename OtherDerived>
-EIGEN_STRONG_INLINE const CwiseBinaryOp<ei_scalar_sum_op<typename ei_traits<Derived>::Scalar>, Derived, OtherDerived>
-MatrixBase<Derived>::operator+(const MatrixBase<OtherDerived> &other) const
-{
- return CwiseBinaryOp<ei_scalar_sum_op<Scalar>, Derived, OtherDerived>(derived(), other.derived());
-}
-
/** replaces \c *this by \c *this + \a other.
*
* \returns a reference to \c *this
@@ -225,111 +193,4 @@ MatrixBase<Derived>::operator+=(const MatrixBase<OtherDerived>& other)
return *this = *this + other;
}
-/** \returns an expression of the Schur product (coefficient wise product) of *this and \a other
- *
- * Example: \include Cwise_product.cpp
- * Output: \verbinclude Cwise_product.out
- *
- * \sa class CwiseBinaryOp, operator/(), square()
- */
-template<typename ExpressionType,template <typename> class StorageBase>
-template<typename OtherDerived>
-EIGEN_STRONG_INLINE const EIGEN_CWISE_PRODUCT_RETURN_TYPE
-Cwise<ExpressionType,StorageBase>::operator*(const AnyMatrixBase<OtherDerived> &other) const
-{
- return EIGEN_CWISE_PRODUCT_RETURN_TYPE(_expression(), other.derived());
-}
-
-/** \returns an expression of the coefficient-wise quotient of *this and \a other
- *
- * Example: \include Cwise_quotient.cpp
- * Output: \verbinclude Cwise_quotient.out
- *
- * \sa class CwiseBinaryOp, operator*(), inverse()
- */
-template<typename ExpressionType,template <typename> class StorageBase>
-template<typename OtherDerived>
-EIGEN_STRONG_INLINE const EIGEN_CWISE_BINOP_RETURN_TYPE(ei_scalar_quotient_op)
-Cwise<ExpressionType,StorageBase>::operator/(const StorageBase<OtherDerived> &other) const
-{
- return EIGEN_CWISE_BINOP_RETURN_TYPE(ei_scalar_quotient_op)(_expression(), other.derived());
-}
-
-/** Replaces this expression by its coefficient-wise product with \a other.
- *
- * Example: \include Cwise_times_equal.cpp
- * Output: \verbinclude Cwise_times_equal.out
- *
- * \sa operator*(), operator/=()
- */
-template<typename ExpressionType,template <typename> class StorageBase>
-template<typename OtherDerived>
-inline ExpressionType& Cwise<ExpressionType,StorageBase>::operator*=(const StorageBase<OtherDerived> &other)
-{
- return m_matrix.const_cast_derived() = *this * other;
-}
-
-/** Replaces this expression by its coefficient-wise quotient by \a other.
- *
- * Example: \include Cwise_slash_equal.cpp
- * Output: \verbinclude Cwise_slash_equal.out
- *
- * \sa operator/(), operator*=()
- */
-template<typename ExpressionType,template <typename> class StorageBase>
-template<typename OtherDerived>
-inline ExpressionType& Cwise<ExpressionType,StorageBase>::operator/=(const StorageBase<OtherDerived> &other)
-{
- return m_matrix.const_cast_derived() = *this / other;
-}
-
-/** \returns an expression of the coefficient-wise min of *this and \a other
- *
- * Example: \include Cwise_min.cpp
- * Output: \verbinclude Cwise_min.out
- *
- * \sa class CwiseBinaryOp
- */
-template<typename ExpressionType,template <typename> class StorageBase>
-template<typename OtherDerived>
-EIGEN_STRONG_INLINE const EIGEN_CWISE_BINOP_RETURN_TYPE(ei_scalar_min_op)
-Cwise<ExpressionType,StorageBase>::min(const StorageBase<OtherDerived> &other) const
-{
- return EIGEN_CWISE_BINOP_RETURN_TYPE(ei_scalar_min_op)(_expression(), other.derived());
-}
-
-/** \returns an expression of the coefficient-wise max of *this and \a other
- *
- * Example: \include Cwise_max.cpp
- * Output: \verbinclude Cwise_max.out
- *
- * \sa class CwiseBinaryOp
- */
-template<typename ExpressionType,template <typename> class StorageBase>
-template<typename OtherDerived>
-EIGEN_STRONG_INLINE const EIGEN_CWISE_BINOP_RETURN_TYPE(ei_scalar_max_op)
-Cwise<ExpressionType,StorageBase>::max(const StorageBase<OtherDerived> &other) const
-{
- return EIGEN_CWISE_BINOP_RETURN_TYPE(ei_scalar_max_op)(_expression(), other.derived());
-}
-
-/** \returns an expression of a custom coefficient-wise operator \a func of *this and \a other
- *
- * The template parameter \a CustomBinaryOp is the type of the functor
- * of the custom operator (see class CwiseBinaryOp for an example)
- *
- * Here is an example illustrating the use of custom functors:
- * \include class_CwiseBinaryOp.cpp
- * Output: \verbinclude class_CwiseBinaryOp.out
- *
- * \sa class CwiseBinaryOp, MatrixBase::operator+, MatrixBase::operator-, Cwise::operator*, Cwise::operator/
- */
-template<typename Derived>
-template<typename CustomBinaryOp, typename OtherDerived>
-EIGEN_STRONG_INLINE const CwiseBinaryOp<CustomBinaryOp, Derived, OtherDerived>
-MatrixBase<Derived>::binaryExpr(const MatrixBase<OtherDerived> &other, const CustomBinaryOp& func) const
-{
- return CwiseBinaryOp<CustomBinaryOp, Derived, OtherDerived>(derived(), other.derived(), func);
-}
-
#endif // EIGEN_CWISE_BINARY_OP_H
diff --git a/Eigen/src/Core/CwiseBinaryOps.h b/Eigen/src/Core/CwiseBinaryOps.h
new file mode 100644
index 000000000..0b7fa2d8b
--- /dev/null
+++ b/Eigen/src/Core/CwiseBinaryOps.h
@@ -0,0 +1,139 @@
+/** \returns an expression of the difference of \c *this and \a other
+ *
+ * \note If you want to substract a given scalar from all coefficients, see Cwise::operator-().
+ *
+ * \sa class CwiseBinaryOp, MatrixBase::operator-=()
+ */
+template<typename OtherDerived>
+EIGEN_STRONG_INLINE const CwiseBinaryOp<ei_scalar_difference_op<typename ei_traits<Derived>::Scalar>,
+ Derived, OtherDerived>
+operator-(const EIGEN_CURRENT_STORAGE_BASE_CLASS<OtherDerived> &other) const
+{
+ return CwiseBinaryOp<ei_scalar_difference_op<Scalar>,
+ Derived, OtherDerived>(derived(), other.derived());
+}
+
+/** \returns an expression of the sum of \c *this and \a other
+ *
+ * \note If you want to add a given scalar to all coefficients, see Cwise::operator+().
+ *
+ * \sa class CwiseBinaryOp, MatrixBase::operator+=()
+ */
+template<typename OtherDerived>
+EIGEN_STRONG_INLINE const CwiseBinaryOp<ei_scalar_sum_op<typename ei_traits<Derived>::Scalar>, Derived, OtherDerived>
+operator+(const EIGEN_CURRENT_STORAGE_BASE_CLASS<OtherDerived> &other) const
+{
+ return CwiseBinaryOp<ei_scalar_sum_op<Scalar>, Derived, OtherDerived>(derived(), other.derived());
+}
+
+/** \returns an expression of a custom coefficient-wise operator \a func of *this and \a other
+ *
+ * The template parameter \a CustomBinaryOp is the type of the functor
+ * of the custom operator (see class CwiseBinaryOp for an example)
+ *
+ * Here is an example illustrating the use of custom functors:
+ * \include class_CwiseBinaryOp.cpp
+ * Output: \verbinclude class_CwiseBinaryOp.out
+ *
+ * \sa class CwiseBinaryOp, MatrixBase::operator+, MatrixBase::operator-, MatrixBase::cwiseProduct
+ */
+template<typename CustomBinaryOp, typename OtherDerived>
+EIGEN_STRONG_INLINE const CwiseBinaryOp<CustomBinaryOp, Derived, OtherDerived>
+binaryExpr(const EIGEN_CURRENT_STORAGE_BASE_CLASS<OtherDerived> &other, const CustomBinaryOp& func = CustomBinaryOp()) const
+{
+ return CwiseBinaryOp<CustomBinaryOp, Derived, OtherDerived>(derived(), other.derived(), func);
+}
+
+/** \returns an expression of the Schur product (coefficient wise product) of *this and \a other
+ *
+ * Example: \include MatrixBase_cwiseProduct.cpp
+ * Output: \verbinclude MatrixBase_cwiseProduct.out
+ *
+ * \sa class CwiseBinaryOp, cwiseAbs2
+ */
+template<typename OtherDerived>
+EIGEN_STRONG_INLINE const CwiseBinaryOp<ei_scalar_product_op<Scalar>, Derived, OtherDerived>
+cwiseProduct(const EIGEN_CURRENT_STORAGE_BASE_CLASS<OtherDerived> &other) const
+{
+ return CwiseBinaryOp<ei_scalar_product_op<Scalar>, Derived, OtherDerived>(derived(), other.derived());
+}
+
+/** \returns an expression of the coefficient-wise == operator of *this and \a other
+ *
+ * \warning this performs an exact comparison, which is generally a bad idea with floating-point types.
+ * In order to check for equality between two vectors or matrices with floating-point coefficients, it is
+ * generally a far better idea to use a fuzzy comparison as provided by MatrixBase::isApprox() and
+ * MatrixBase::isMuchSmallerThan().
+ *
+ * Example: \include MatrixBase_cwiseEqual.cpp
+ * Output: \verbinclude MatrixBase_cwiseEqual.out
+ *
+ * \sa MatrixBase::cwiseNotEqual(), MatrixBase::isApprox(), MatrixBase::isMuchSmallerThan()
+ */
+template<typename OtherDerived>
+inline const CwiseBinaryOp<std::equal_to<Scalar>, Derived, OtherDerived>
+cwiseEqual(const EIGEN_CURRENT_STORAGE_BASE_CLASS<OtherDerived> &other) const
+{
+ return CwiseBinaryOp<std::equal_to<Scalar>, Derived, OtherDerived>(derived(), other.derived());
+}
+
+/** \returns an expression of the coefficient-wise != operator of *this and \a other
+ *
+ * \warning this performs an exact comparison, which is generally a bad idea with floating-point types.
+ * In order to check for equality between two vectors or matrices with floating-point coefficients, it is
+ * generally a far better idea to use a fuzzy comparison as provided by MatrixBase::isApprox() and
+ * MatrixBase::isMuchSmallerThan().
+ *
+ * Example: \include MatrixBase_cwiseNotEqual.cpp
+ * Output: \verbinclude MatrixBase_cwiseNotEqual.out
+ *
+ * \sa MatrixBase::cwiseEqual(), MatrixBase::isApprox(), MatrixBase::isMuchSmallerThan()
+ */
+template<typename OtherDerived>
+inline const CwiseBinaryOp<std::not_equal_to<Scalar>, Derived, OtherDerived>
+cwiseNotEqual(const EIGEN_CURRENT_STORAGE_BASE_CLASS<OtherDerived> &other) const
+{
+ return CwiseBinaryOp<std::not_equal_to<Scalar>, Derived, OtherDerived>(derived(), other.derived());
+}
+
+/** \returns an expression of the coefficient-wise min of *this and \a other
+ *
+ * Example: \include MatrixBase_cwiseMin.cpp
+ * Output: \verbinclude MatrixBase_cwiseMin.out
+ *
+ * \sa class CwiseBinaryOp, max()
+ */
+template<typename OtherDerived>
+EIGEN_STRONG_INLINE const CwiseBinaryOp<ei_scalar_min_op<Scalar>, Derived, OtherDerived>
+cwiseMin(const EIGEN_CURRENT_STORAGE_BASE_CLASS<OtherDerived> &other) const
+{
+ return CwiseBinaryOp<ei_scalar_min_op<Scalar>, Derived, OtherDerived>(derived(), other.derived());
+}
+
+/** \returns an expression of the coefficient-wise max of *this and \a other
+ *
+ * Example: \include MatrixBase_cwiseMax.cpp
+ * Output: \verbinclude MatrixBase_cwiseMax.out
+ *
+ * \sa class CwiseBinaryOp, min()
+ */
+template<typename OtherDerived>
+EIGEN_STRONG_INLINE const CwiseBinaryOp<ei_scalar_max_op<Scalar>, Derived, OtherDerived>
+cwiseMax(const EIGEN_CURRENT_STORAGE_BASE_CLASS<OtherDerived> &other) const
+{
+ return CwiseBinaryOp<ei_scalar_max_op<Scalar>, Derived, OtherDerived>(derived(), other.derived());
+}
+
+/** \returns an expression of the coefficient-wise quotient of *this and \a other
+ *
+ * Example: \include MatrixBase_cwiseQuotient.cpp
+ * Output: \verbinclude MatrixBase_cwiseQuotient.out
+ *
+ * \sa class CwiseBinaryOp, cwiseProduct(), cwiseInverse()
+ */
+template<typename OtherDerived>
+EIGEN_STRONG_INLINE const CwiseBinaryOp<ei_scalar_quotient_op<Scalar>, Derived, OtherDerived>
+cwiseQuotient(const EIGEN_CURRENT_STORAGE_BASE_CLASS<OtherDerived> &other) const
+{
+ return CwiseBinaryOp<ei_scalar_quotient_op<Scalar>, Derived, OtherDerived>(derived(), other.derived());
+}
diff --git a/Eigen/src/Core/CwiseUnaryOp.h b/Eigen/src/Core/CwiseUnaryOp.h
index 0ba7e1366..55f965b4e 100644
--- a/Eigen/src/Core/CwiseUnaryOp.h
+++ b/Eigen/src/Core/CwiseUnaryOp.h
@@ -128,60 +128,4 @@ class CwiseUnaryOpImpl<UnaryOp,MatrixType,Dense> : public MatrixBase<CwiseUnaryO
}
};
-/** \returns an expression of the coefficient-wise absolute value of \c *this
- *
- * Example: \include Cwise_abs.cpp
- * Output: \verbinclude Cwise_abs.out
- *
- * \sa abs2()
- */
-template<typename ExpressionType,template <typename> class StorageBase>
-EIGEN_STRONG_INLINE const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_abs_op)
-Cwise<ExpressionType,StorageBase>::abs() const
-{
- return _expression();
-}
-
-/** \returns an expression of the coefficient-wise squared absolute value of \c *this
- *
- * Example: \include Cwise_abs2.cpp
- * Output: \verbinclude Cwise_abs2.out
- *
- * \sa abs(), square()
- */
-template<typename ExpressionType,template <typename> class StorageBase>
-EIGEN_STRONG_INLINE const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_abs2_op)
-Cwise<ExpressionType,StorageBase>::abs2() const
-{
- return _expression();
-}
-
-/** \returns an expression of the coefficient-wise exponential of *this.
- *
- * Example: \include Cwise_exp.cpp
- * Output: \verbinclude Cwise_exp.out
- *
- * \sa pow(), log(), sin(), cos()
- */
-template<typename ExpressionType,template <typename> class StorageBase>
-inline const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_exp_op)
-Cwise<ExpressionType,StorageBase>::exp() const
-{
- return _expression();
-}
-
-/** \returns an expression of the coefficient-wise logarithm of *this.
- *
- * Example: \include Cwise_log.cpp
- * Output: \verbinclude Cwise_log.out
- *
- * \sa exp()
- */
-template<typename ExpressionType,template <typename> class StorageBase>
-inline const EIGEN_CWISE_UNOP_RETURN_TYPE(ei_scalar_log_op)
-Cwise<ExpressionType,StorageBase>::log() const
-{
- return _expression();
-}
-
#endif // EIGEN_CWISE_UNARY_OP_H
diff --git a/Eigen/src/Core/CwiseUnaryOps.h b/Eigen/src/Core/CwiseUnaryOps.h
index cd58f1d43..39fd479b5 100644
--- a/Eigen/src/Core/CwiseUnaryOps.h
+++ b/Eigen/src/Core/CwiseUnaryOps.h
@@ -156,26 +156,58 @@ real() { return derived(); }
EIGEN_STRONG_INLINE NonConstImagReturnType
imag() { return derived(); }
-/** \returns a Cwise wrapper of *this providing additional coefficient-wise operations
+/** \returns an expression of the coefficient-wise absolute value of \c *this
*
- * Example: \include MatrixBase_cwise_const.cpp
- * Output: \verbinclude MatrixBase_cwise_const.out
+ * Example: \include MatrixBase_cwiseAbs.cpp
+ * Output: \verbinclude MatrixBase_cwiseAbs.out
*
- * \sa class Cwise, cwise()
+ * \sa cwiseAbs2()
*/
-inline const Cwise<Derived,EIGEN_CURRENT_STORAGE_BASE_CLASS> cwise() const
-{
- return derived();
-}
+EIGEN_STRONG_INLINE const CwiseUnaryOp<ei_scalar_abs_op<Scalar>,Derived>
+cwiseAbs() const { return derived(); }
+
+/** \returns an expression of the coefficient-wise squared absolute value of \c *this
+ *
+ * Example: \include MatrixBase_cwiseAbs2.cpp
+ * Output: \verbinclude MatrixBase_cwiseAbs2.out
+ *
+ * \sa cwiseAbs()
+ */
+EIGEN_STRONG_INLINE const CwiseUnaryOp<ei_scalar_abs2_op<Scalar>,Derived>
+cwiseAbs2() const { return derived(); }
+
+/** \returns an expression of the coefficient-wise square root of *this.
+ *
+ * Example: \include MatrixBase_cwiseSqrt.cpp
+ * Output: \verbinclude MatrixBase_cwiseSqrt.out
+ *
+ * \sa cwisePow(), cwiseSquare()
+ */
+inline const CwiseUnaryOp<ei_scalar_sqrt_op<Scalar>,Derived>
+cwiseSqrt() const { return derived(); }
+
+/** \returns an expression of the coefficient-wise inverse of *this.
+ *
+ * Example: \include MatrixBase_cwiseInverse.cpp
+ * Output: \verbinclude MatrixBase_cwiseInverse.out
+ *
+ * \sa cwiseProduct()
+ */
+inline const CwiseUnaryOp<ei_scalar_inverse_op<Scalar>,Derived>
+cwiseInverse() const { return derived(); }
-/** \returns a Cwise wrapper of *this providing additional coefficient-wise operations
+/** \returns an expression of the coefficient-wise == operator of \c *this and a scalar \a s
*
- * Example: \include MatrixBase_cwise.cpp
- * Output: \verbinclude MatrixBase_cwise.out
+ * \warning this performs an exact comparison, which is generally a bad idea with floating-point types.
+ * In order to check for equality between two vectors or matrices with floating-point coefficients, it is
+ * generally a far better idea to use a fuzzy comparison as provided by MatrixBase::isApprox() and
+ * MatrixBase::isMuchSmallerThan().
*
- * \sa class Cwise, cwise() const
+ * \sa cwiseEqual(const MatrixBase<OtherDerived> &) const
*/
-inline Cwise<Derived,EIGEN_CURRENT_STORAGE_BASE_CLASS> cwise()
+inline const CwiseUnaryOp<std::binder1st<std::equal_to<Scalar> >,Derived>
+cwiseEqual(Scalar s) const
{
- return derived();
+ return CwiseUnaryOp<std::binder1st<std::equal_to<Scalar> >,Derived>
+ (derived(), std::bind1st(std::equal_to<Scalar>(), s));
}
diff --git a/Eigen/src/Core/DiagonalMatrix.h b/Eigen/src/Core/DiagonalMatrix.h
index 1dec82229..6f93737ff 100644
--- a/Eigen/src/Core/DiagonalMatrix.h
+++ b/Eigen/src/Core/DiagonalMatrix.h
@@ -68,6 +68,12 @@ class DiagonalBase : public AnyMatrixBase<Derived>
template<typename MatrixDerived>
const DiagonalProduct<MatrixDerived, Derived, OnTheLeft>
operator*(const MatrixBase<MatrixDerived> &matrix) const;
+
+ inline const DiagonalWrapper<NestByValue<CwiseUnaryOp<ei_scalar_inverse_op<Scalar>, DiagonalVectorType> > >
+ inverse() const
+ {
+ return diagonal().cwiseInverse().nestByValue();
+ }
};
template<typename Derived>
diff --git a/Eigen/src/Core/Dot.h b/Eigen/src/Core/Dot.h
index 631124f2b..4a164a99e 100644
--- a/Eigen/src/Core/Dot.h
+++ b/Eigen/src/Core/Dot.h
@@ -279,7 +279,7 @@ MatrixBase<Derived>::dot(const MatrixBase<OtherDerived>& other) const
template<typename Derived>
inline typename NumTraits<typename ei_traits<Derived>::Scalar>::Real MatrixBase<Derived>::squaredNorm() const
{
- return ei_real((*this).cwise().abs2().sum());
+ return ei_real((*this).cwiseAbs2().sum());
}
/** \returns the \em l2 norm of *this, i.e., for vectors, the square root of the dot product of *this with itself.
diff --git a/Eigen/src/Core/Flagged.h b/Eigen/src/Core/Flagged.h
deleted file mode 100644
index 754eaf6c5..000000000
--- a/Eigen/src/Core/Flagged.h
+++ /dev/null
@@ -1,149 +0,0 @@
-// This file is part of Eigen, a lightweight C++ template library
-// for linear algebra.
-//
-// Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
-//
-// Eigen is free software; you can redistribute it and/or
-// modify it under the terms of the GNU Lesser General Public
-// License as published by the Free Software Foundation; either
-// version 3 of the License, or (at your option) any later version.
-//
-// Alternatively, you can redistribute it and/or
-// modify it under the terms of the GNU General Public License as
-// published by the Free Software Foundation; either version 2 of
-// the License, or (at your option) any later version.
-//
-// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
-// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
-// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
-// GNU General Public License for more details.
-//
-// You should have received a copy of the GNU Lesser General Public
-// License and a copy of the GNU General Public License along with
-// Eigen. If not, see <http://www.gnu.org/licenses/>.
-
-#ifndef EIGEN_FLAGGED_H
-#define EIGEN_FLAGGED_H
-
-/** \deprecated it is only used by lazy() which is deprecated
- *
- * \class Flagged
- *
- * \brief Expression with modified flags
- *
- * \param ExpressionType the type of the object of which we are modifying the flags
- * \param Added the flags added to the expression
- * \param Removed the flags removed from the expression (has priority over Added).
- *
- * This class represents an expression whose flags have been modified.
- * It is the return type of MatrixBase::flagged()
- * and most of the time this is the only way it is used.
- *
- * \sa MatrixBase::flagged()
- */
-template<typename ExpressionType, unsigned int Added, unsigned int Removed>
-struct ei_traits<Flagged<ExpressionType, Added, Removed> > : ei_traits<ExpressionType>
-{
- enum { Flags = (ExpressionType::Flags | Added) & ~Removed };
-};
-
-template<typename ExpressionType, unsigned int Added, unsigned int Removed> class Flagged
- : public MatrixBase<Flagged<ExpressionType, Added, Removed> >
-{
- public:
-
- EIGEN_GENERIC_PUBLIC_INTERFACE(Flagged)
- typedef typename ei_meta_if<ei_must_nest_by_value<ExpressionType>::ret,
- ExpressionType, const ExpressionType&>::ret ExpressionTypeNested;
- typedef typename ExpressionType::InnerIterator InnerIterator;
-
- inline Flagged(const ExpressionType& matrix) : m_matrix(matrix) {}
-
- inline int rows() const { return m_matrix.rows(); }
- inline int cols() const { return m_matrix.cols(); }
- inline int stride() const { return m_matrix.stride(); }
-
- inline const Scalar coeff(int row, int col) const
- {
- return m_matrix.coeff(row, col);
- }
-
- inline Scalar& coeffRef(int row, int col)
- {
- return m_matrix.const_cast_derived().coeffRef(row, col);
- }
-
- inline const Scalar coeff(int index) const
- {
- return m_matrix.coeff(index);
- }
-
- inline Scalar& coeffRef(int index)
- {
- return m_matrix.const_cast_derived().coeffRef(index);
- }
-
- template<int LoadMode>
- inline const PacketScalar packet(int row, int col) const
- {
- return m_matrix.template packet<LoadMode>(row, col);
- }
-
- template<int LoadMode>
- inline void writePacket(int row, int col, const PacketScalar& x)
- {
- m_matrix.const_cast_derived().template writePacket<LoadMode>(row, col, x);
- }
-
- template<int LoadMode>
- inline const PacketScalar packet(int index) const
- {
- return m_matrix.template packet<LoadMode>(index);
- }
-
- template<int LoadMode>
- inline void writePacket(int index, const PacketScalar& x)
- {
- m_matrix.const_cast_derived().template writePacket<LoadMode>(index, x);
- }
-
- const ExpressionType& _expression() const { return m_matrix; }
-
- protected:
- ExpressionTypeNested m_matrix;
-};
-
-/** \deprecated it is only used by lazy() which is deprecated
- *
- * \returns an expression of *this with added flags
- *
- * Example: \include MatrixBase_marked.cpp
- * Output: \verbinclude MatrixBase_marked.out
- *
- * \sa class Flagged, extract(), part()
- */
-template<typename Derived>
-template<unsigned int Added>
-inline const Flagged<Derived, Added, 0>
-MatrixBase<Derived>::marked() const
-{
- return derived();
-}
-
-/** \deprecated use MatrixBase::noalias()
- *
- * \returns an expression of *this with the EvalBeforeAssigningBit flag removed.
- *
- * Example: \include MatrixBase_lazy.cpp
- * Output: \verbinclude MatrixBase_lazy.out
- *
- * \sa class Flagged, marked()
- */
-template<typename Derived>
-inline const Flagged<Derived, 0, EvalBeforeAssigningBit>
-MatrixBase<Derived>::lazy() const
-{
- return derived();
-}
-
-#endif // EIGEN_FLAGGED_H
diff --git a/Eigen/src/Core/Fuzzy.h b/Eigen/src/Core/Fuzzy.h
index e10446398..13fcfcdab 100644
--- a/Eigen/src/Core/Fuzzy.h
+++ b/Eigen/src/Core/Fuzzy.h
@@ -54,7 +54,7 @@ bool MatrixBase<Derived>::isApprox(
{
const typename ei_nested<Derived,2>::type nested(derived());
const typename ei_nested<OtherDerived,2>::type otherNested(other.derived());
- return (nested - otherNested).cwise().abs2().sum() <= prec * prec * std::min(nested.cwise().abs2().sum(), otherNested.cwise().abs2().sum());
+ return (nested - otherNested).cwiseAbs2().sum() <= prec * prec * std::min(nested.cwiseAbs2().sum(), otherNested.cwiseAbs2().sum());
}
/** \returns \c true if the norm of \c *this is much smaller than \a other,
@@ -76,7 +76,7 @@ bool MatrixBase<Derived>::isMuchSmallerThan(
RealScalar prec
) const
{
- return cwise().abs2().sum() <= prec * prec * other * other;
+ return cwiseAbs2().sum() <= prec * prec * other * other;
}
/** \returns \c true if the norm of \c *this is much smaller than the norm of \a other,
@@ -96,7 +96,7 @@ bool MatrixBase<Derived>::isMuchSmallerThan(
RealScalar prec
) const
{
- return this->cwise().abs2().sum() <= prec * prec * other.cwise().abs2().sum();
+ return cwiseAbs2().sum() <= prec * prec * other.cwiseAbs2().sum();
}
#else
diff --git a/Eigen/src/Core/MatrixBase.h b/Eigen/src/Core/MatrixBase.h
index 67e270af7..d653b25c5 100644
--- a/Eigen/src/Core/MatrixBase.h
+++ b/Eigen/src/Core/MatrixBase.h
@@ -247,6 +247,7 @@ template<typename Derived> class MatrixBase
#define EIGEN_CURRENT_STORAGE_BASE_CLASS Eigen::MatrixBase
#include "CwiseUnaryOps.h"
+ #include "CwiseBinaryOps.h"
#undef EIGEN_CURRENT_STORAGE_BASE_CLASS
/** Copies \a other into *this. \returns a reference to *this. */
@@ -275,12 +276,6 @@ template<typename Derived> class MatrixBase
template<typename OtherDerived>
Derived& lazyAssign(const MatrixBase<OtherDerived>& other);
- /** \deprecated because .lazy() is deprecated
- * Overloaded for cache friendly product evaluation */
- template<typename OtherDerived>
- Derived& lazyAssign(const Flagged<OtherDerived, 0, EvalBeforeAssigningBit>& other)
- { return lazyAssign(other._expression()); }
-
template<typename ProductDerived, typename Lhs, typename Rhs>
Derived& lazyAssign(const ProductBase<ProductDerived, Lhs,Rhs>& other);
@@ -342,15 +337,6 @@ template<typename Derived> class MatrixBase
Scalar& z();
Scalar& w();
-
- template<typename OtherDerived>
- const CwiseBinaryOp<ei_scalar_sum_op<typename ei_traits<Derived>::Scalar>, Derived, OtherDerived>
- operator+(const MatrixBase<OtherDerived> &other) const;
-
- template<typename OtherDerived>
- const CwiseBinaryOp<ei_scalar_difference_op<typename ei_traits<Derived>::Scalar>, Derived, OtherDerived>
- operator-(const MatrixBase<OtherDerived> &other) const;
-
template<typename OtherDerived>
Derived& operator+=(const MatrixBase<OtherDerived>& other);
template<typename OtherDerived>
@@ -374,14 +360,6 @@ template<typename Derived> class MatrixBase
operator*(const DiagonalBase<DiagonalDerived> &diagonal) const;
template<typename OtherDerived>
- typename ei_plain_matrix_type_column_major<OtherDerived>::type
- solveTriangular(const MatrixBase<OtherDerived>& other) const;
-
- template<typename OtherDerived>
- void solveTriangularInPlace(const MatrixBase<OtherDerived>& other) const;
-
-
- template<typename OtherDerived>
Scalar dot(const MatrixBase<OtherDerived>& other) const;
RealScalar squaredNorm() const;
RealScalar norm() const;
@@ -542,11 +520,11 @@ template<typename Derived> class MatrixBase
template<typename OtherDerived>
inline bool operator==(const MatrixBase<OtherDerived>& other) const
- { return (cwise() == other).all(); }
+ { return cwiseEqual(other).all(); }
template<typename OtherDerived>
inline bool operator!=(const MatrixBase<OtherDerived>& other) const
- { return (cwise() != other).any(); }
+ { return cwiseNotEqual(other).all(); }
/** \returns the matrix or vector obtained by evaluating this expression.
@@ -560,10 +538,6 @@ template<typename Derived> class MatrixBase
template<typename OtherDerived>
void swap(MatrixBase<OtherDerived> EIGEN_REF_TO_TEMPORARY other);
- template<unsigned int Added>
- const Flagged<Derived, Added, 0> marked() const;
- const Flagged<Derived, 0, EvalBeforeAssigningBit> lazy() const;
-
NoAlias<Derived,Eigen::MatrixBase > noalias();
/** \returns number of elements to skip to pass from one row (resp. column) to another
@@ -575,12 +549,6 @@ template<typename Derived> class MatrixBase
inline const NestByValue<Derived> nestByValue() const;
-
- template<typename CustomBinaryOp, typename OtherDerived>
- const CwiseBinaryOp<CustomBinaryOp, Derived, OtherDerived>
- binaryExpr(const MatrixBase<OtherDerived> &other, const CustomBinaryOp& func = CustomBinaryOp()) const;
-
-
Scalar sum() const;
Scalar mean() const;
Scalar trace() const;
@@ -736,6 +704,33 @@ template<typename Derived> class MatrixBase
#include EIGEN_MATRIXBASE_PLUGIN
#endif
+#ifdef EIGEN2_SUPPORT
+ /** \deprecated because .lazy() is deprecated
+ * Overloaded for cache friendly product evaluation */
+ template<typename OtherDerived>
+ Derived& lazyAssign(const Flagged<OtherDerived, 0, EvalBeforeAssigningBit>& other)
+ { return lazyAssign(other._expression()); }
+
+ template<unsigned int Added>
+ const Flagged<Derived, Added, 0> marked() const;
+ const Flagged<Derived, 0, EvalBeforeAssigningBit> lazy() const;
+
+ inline const Cwise<Derived> cwise() const;
+ inline Cwise<Derived> cwise();
+
+ // a workaround waiting the Array class
+ inline const Cwise<Derived> array() const { return cwise(); }
+ // a workaround waiting the Array class
+ inline Cwise<Derived> array() { return cwise(); }
+
+ template<typename OtherDerived>
+ typename ei_plain_matrix_type_column_major<OtherDerived>::type
+ solveTriangular(const MatrixBase<OtherDerived>& other) const;
+
+ template<typename OtherDerived>
+ void solveTriangularInPlace(const MatrixBase<OtherDerived>& other) const;
+#endif
+
protected:
/** Default constructor. Do nothing. */
MatrixBase()
diff --git a/Eigen/src/Core/Product.h b/Eigen/src/Core/Product.h
index cc751650d..7db35eaad 100644
--- a/Eigen/src/Core/Product.h
+++ b/Eigen/src/Core/Product.h
@@ -170,7 +170,7 @@ class GeneralProduct<Lhs, Rhs, InnerProduct>
EIGEN_STRONG_INLINE Scalar value() const
{
- return (m_lhs.transpose().cwise()*m_rhs).sum();
+ return (m_lhs.transpose().cwiseProduct(m_rhs)).sum();
}
template<typename Dest> void scaleAndAddTo(Dest& dst, Scalar alpha) const
@@ -403,7 +403,7 @@ template<> struct ei_gemv_selector<OnTheRight,RowMajor,false>
// TODO makes sure rhs is sequentially stored in memory, otherwise use a temp
const int rows = prod.rows();
for(int i=0; i<rows; ++i)
- dest.coeffRef(i) += alpha * (prod.lhs().row(i).cwise() * prod.rhs().transpose()).sum();
+ dest.coeffRef(i) += alpha * (prod.lhs().row(i).cwiseProduct(prod.rhs().transpose())).sum();
}
};
@@ -431,7 +431,7 @@ MatrixBase<Derived>::operator*(const MatrixBase<OtherDerived> &other) const
};
// note to the lost user:
// * for a dot product use: v1.dot(v2)
- // * for a coeff-wise product use: v1.cwise()*v2
+ // * for a coeff-wise product use: v1.cwiseProduct(v2)
EIGEN_STATIC_ASSERT(ProductIsValid || !(AreVectors && SameSizes),
INVALID_VECTOR_VECTOR_PRODUCT__IF_YOU_WANTED_A_DOT_OR_COEFF_WISE_PRODUCT_YOU_MUST_USE_THE_EXPLICIT_FUNCTIONS)
EIGEN_STATIC_ASSERT(ProductIsValid || !(SameSizes && !AreVectors),
diff --git a/Eigen/src/Core/SolveTriangular.h b/Eigen/src/Core/SolveTriangular.h
index c7f0cd227..6e9dd30ec 100644
--- a/Eigen/src/Core/SolveTriangular.h
+++ b/Eigen/src/Core/SolveTriangular.h
@@ -78,8 +78,7 @@ struct ei_triangular_solver_selector<Lhs,Rhs,OnTheLeft,Mode,NoUnrolling,RowMajor
int i = IsLowerTriangular ? pi+k : pi-k-1;
int s = IsLowerTriangular ? pi : i+1;
if (k>0)
- other.coeffRef(i) -= ((lhs.row(i).segment(s,k).transpose())
- .cwise()*(other.segment(s,k))).sum();
+ other.coeffRef(i) -= (lhs.row(i).segment(s,k).transpose().cwiseProduct(other.segment(s,k))).sum();
if(!(Mode & UnitDiagBit))
other.coeffRef(i) /= lhs.coeff(i,i);
@@ -180,8 +179,7 @@ struct ei_triangular_solver_unroller<Lhs,Rhs,Mode,Index,Size,false> {
static void run(const Lhs& lhs, Rhs& rhs)
{
if (Index>0)
- rhs.coeffRef(I) -= ((lhs.row(I).template segment<Index>(S).transpose())
- .cwise()*(rhs.template segment<Index>(S))).sum();
+ rhs.coeffRef(I) -= ((lhs.row(I).template segment<Index>(S).transpose()).cwiseProduct(rhs.template segment<Index>(S))).sum();
if(!(Mode & UnitDiagBit))
rhs.coeffRef(I) /= lhs.coeff(I,I);
diff --git a/Eigen/src/Core/StableNorm.h b/Eigen/src/Core/StableNorm.h
index f2d1e7240..b08f4a1ae 100644
--- a/Eigen/src/Core/StableNorm.h
+++ b/Eigen/src/Core/StableNorm.h
@@ -28,7 +28,7 @@
template<typename ExpressionType, typename Scalar>
inline void ei_stable_norm_kernel(const ExpressionType& bl, Scalar& ssq, Scalar& scale, Scalar& invScale)
{
- Scalar max = bl.cwise().abs().maxCoeff();
+ Scalar max = bl.cwiseAbs().maxCoeff();
if (max>scale)
{
ssq = ssq * ei_abs2(scale/max);
@@ -182,7 +182,7 @@ template<typename Derived>
inline typename NumTraits<typename ei_traits<Derived>::Scalar>::Real
MatrixBase<Derived>::hypotNorm() const
{
- return this->cwise().abs().redux(ei_scalar_hypot_op<RealScalar>());
+ return this->cwiseAbs().redux(ei_scalar_hypot_op<RealScalar>());
}
#endif // EIGEN_STABLENORM_H
diff --git a/Eigen/src/Core/util/ForwardDeclarations.h b/Eigen/src/Core/util/ForwardDeclarations.h
index 2fedbbc07..541b5dd9f 100644
--- a/Eigen/src/Core/util/ForwardDeclarations.h
+++ b/Eigen/src/Core/util/ForwardDeclarations.h
@@ -60,7 +60,6 @@ template<typename MatrixType, int PacketAccess = AsRequested> class Map;
template<typename Derived> class TriangularBase;
template<typename MatrixType, unsigned int Mode> class TriangularView;
template<typename MatrixType, unsigned int Mode> class SelfAdjointView;
-template<typename ExpressionType, template <typename> class StorageBase> class Cwise;
template<typename ExpressionType> class WithFormat;
template<typename MatrixType> struct CommaInitializer;
template<typename Derived> class ReturnByValue;
@@ -146,4 +145,8 @@ template<typename Scalar,int Dim> class Translation;
template<typename Scalar> class UniformScaling;
template<typename MatrixType,int Direction> class Homogeneous;
+#ifdef EIGEN2_SUPPORT
+template<typename ExpressionType> class Cwise;
+#endif
+
#endif // EIGEN_FORWARDDECLARATIONS_H