aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2008-06-20 12:38:03 +0000
committerGravatar Gael Guennebaud <g.gael@free.fr>2008-06-20 12:38:03 +0000
commit54238961d65ba96f1352feca1b3ae20d09362ce7 (patch)
tree2b8cea91f52bd2412ffa7978972f467e009f0da5
parente735692e3709d7b42e1c1c7dfb1055081b49b0dc (diff)
* added a pseudo expression Array giving access to:
- matrix-scalar addition/subtraction operators, e.g.: m.array() += 0.5; - matrix/matrix comparison operators, e.g.: if (m1.array() < m2.array()) {} * fix compilation issues with Transform and gcc < 4.1
-rw-r--r--Eigen/Array1
-rw-r--r--Eigen/src/Array/Array.h152
-rw-r--r--Eigen/src/Array/Functors.h42
-rw-r--r--Eigen/src/Core/MatrixBase.h3
-rw-r--r--Eigen/src/Core/util/ForwardDeclarations.h2
-rw-r--r--Eigen/src/Geometry/Transform.h110
-rw-r--r--test/CMakeLists.txt1
-rw-r--r--test/linearstructure.cpp1
8 files changed, 255 insertions, 57 deletions
diff --git a/Eigen/Array b/Eigen/Array
index 51c3abe31..33ae6d404 100644
--- a/Eigen/Array
+++ b/Eigen/Array
@@ -10,6 +10,7 @@ namespace Eigen {
#include "src/Array/AllAndAny.h"
#include "src/Array/PartialRedux.h"
#include "src/Array/Random.h"
+#include "src/Array/Array.h"
} // namespace Eigen
diff --git a/Eigen/src/Array/Array.h b/Eigen/src/Array/Array.h
new file mode 100644
index 000000000..5fd074e76
--- /dev/null
+++ b/Eigen/src/Array/Array.h
@@ -0,0 +1,152 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra. Eigen itself is part of the KDE project.
+//
+// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
+//
+// 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_ARRAY_H
+#define EIGEN_ARRAY_H
+
+/** \class Array
+ *
+ * \brief Pseudo expression offering additional features to an expression
+ *
+ * \param ExpressionType the type of the object of which we want array related features
+ *
+ * This class represents an expression with additional, array related, features.
+ * It is the return type of MatrixBase::array()
+ * and most of the time this is the only way it is used.
+ *
+ * \sa MatrixBase::array()
+ */
+template<typename ExpressionType> class Array
+{
+ 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 NestByValue<typename ExpressionType::ConstantReturnType> ConstantReturnType;
+ typedef CwiseUnaryOp<ei_scalar_add_op<Scalar>, ExpressionType> ScalarAddReturnType;
+
+ inline Array(const ExpressionType& matrix) : m_matrix(matrix) {}
+
+ /** \internal */
+ inline const ExpressionType& _expression() const { return m_matrix; }
+
+ const ScalarAddReturnType
+ operator+(const Scalar& scalar) const;
+
+ /** \relates Array */
+ friend const ScalarAddReturnType
+ operator+(const Scalar& scalar, const Array& mat)
+ { return mat + scalar; }
+
+ ExpressionType& operator+=(const Scalar& scalar);
+
+ const ScalarAddReturnType
+ operator-(const Scalar& scalar) const;
+
+ ExpressionType& operator-=(const Scalar& scalar);
+
+ /** \returns true if each coeff of \c *this is less than its respective coeff of \a other */
+ template<typename OtherDerived> bool operator<(const Array<OtherDerived>& other) const
+ { return m_matrix.cwiseLessThan(other._expression()).all(); }
+
+ /** \returns true if each coeff of \c *this is less or equal to its respective coeff of \a other */
+ template<typename OtherDerived> bool operator<=(const Array<OtherDerived>& other) const
+ { return m_matrix.cwiseLessEqual(other._expression()).all(); }
+
+ /** \returns true if each coeff of \c *this is greater to its respective coeff of \a other */
+ template<typename OtherDerived>
+ bool operator>(const Array<OtherDerived>& other) const
+ { return m_matrix.cwiseGreaterThan(other._expression()).all(); }
+
+ /** \returns true if each coeff of \c *this is greater or equal to its respective coeff of \a other */
+ template<typename OtherDerived>
+ bool operator>=(const Array<OtherDerived>& other) const
+ { return m_matrix.cwiseGreaterEqual(other._expression()).all(); }
+
+ protected:
+ ExpressionTypeNested m_matrix;
+};
+
+/** \returns an expression of \c *this with each coeff incremented by the constant \a scalar */
+template<typename ExpressionType>
+const typename Array<ExpressionType>::ScalarAddReturnType
+Array<ExpressionType>::operator+(const Scalar& scalar) const
+{
+ return CwiseUnaryOp<ei_scalar_add_op<Scalar>, ExpressionType>(m_matrix, ei_scalar_add_op<Scalar>(scalar));
+}
+
+/** \see operator+ */
+template<typename ExpressionType>
+ExpressionType& Array<ExpressionType>::operator+=(const Scalar& scalar)
+{
+ m_matrix.const_cast_derived() = *this + scalar;
+ return m_matrix.const_cast_derived();
+}
+
+/** \returns an expression of \c *this with each coeff decremented by the constant \a scalar */
+template<typename ExpressionType>
+const typename Array<ExpressionType>::ScalarAddReturnType
+Array<ExpressionType>::operator-(const Scalar& scalar) const
+{
+ return *this + (-scalar);
+}
+
+/** \see operator- */
+template<typename ExpressionType>
+ExpressionType& Array<ExpressionType>::operator-=(const Scalar& scalar)
+{
+ m_matrix.const_cast_derived() = *this - scalar;
+ return m_matrix.const_cast_derived();
+}
+
+/** \array_module
+ *
+ * \returns an Array expression of *this providing additional,
+ * array related, features.
+ *
+ * \sa class Array
+ */
+template<typename Derived>
+inline const Array<Derived>
+MatrixBase<Derived>::array() const
+{
+ return derived();
+}
+
+/** \array_module
+ *
+ * \returns an Array expression of *this providing additional,
+ * array related, features.
+ *
+ * \sa class Array
+ */
+template<typename Derived>
+inline Array<Derived>
+MatrixBase<Derived>::array()
+{
+ return derived();
+}
+
+#endif // EIGEN_FLAGGED_H
diff --git a/Eigen/src/Array/Functors.h b/Eigen/src/Array/Functors.h
index 6a0a25a3a..03ae98090 100644
--- a/Eigen/src/Array/Functors.h
+++ b/Eigen/src/Array/Functors.h
@@ -26,9 +26,37 @@
#define EIGEN_ARRAY_FUNCTORS_H
/** \internal
+ * \array_module
+ *
+ * \brief Template functor to add a scalar to a fixed other one
+ *
+ * \sa class CwiseUnaryOp, Array::operator+
+ */
+template<typename Scalar, bool PacketAccess = (int(ei_packet_traits<Scalar>::size)>1?true:false) > struct ei_scalar_add_op;
+
+template<typename Scalar>
+struct ei_scalar_add_op<Scalar,true> {
+ typedef typename ei_packet_traits<Scalar>::type PacketScalar;
+ inline ei_scalar_add_op(const Scalar& other) : m_other(ei_pset1(other)) { }
+ inline Scalar operator() (const Scalar& a) const { return a + ei_pfirst(m_other); }
+ inline const PacketScalar packetOp(const PacketScalar& a) const
+ { return ei_padd(a, m_other); }
+ const PacketScalar m_other;
+};
+template<typename Scalar>
+struct ei_scalar_add_op<Scalar,false> {
+ inline ei_scalar_add_op(const Scalar& other) : m_other(other) { }
+ inline Scalar operator() (const Scalar& a) const { return a + m_other; }
+ const Scalar m_other;
+};
+template<typename Scalar>
+struct ei_functor_traits<ei_scalar_add_op<Scalar> >
+{ enum { Cost = NumTraits<Scalar>::AddCost, PacketAccess = ei_packet_traits<Scalar>::size>1 }; };
+
+/** \internal
*
* \array_module
- *
+ *
* \brief Template functor to compute the square root of a scalar
*
* \sa class CwiseUnaryOp, MatrixBase::cwiseSqrt()
@@ -43,7 +71,7 @@ struct ei_functor_traits<ei_scalar_sqrt_op<Scalar> >
/** \internal
*
* \array_module
- *
+ *
* \brief Template functor to compute the exponential of a scalar
*
* \sa class CwiseUnaryOp, MatrixBase::cwiseExp()
@@ -58,7 +86,7 @@ struct ei_functor_traits<ei_scalar_exp_op<Scalar> >
/** \internal
*
* \array_module
- *
+ *
* \brief Template functor to compute the logarithm of a scalar
*
* \sa class CwiseUnaryOp, MatrixBase::cwiseLog()
@@ -73,7 +101,7 @@ struct ei_functor_traits<ei_scalar_log_op<Scalar> >
/** \internal
*
* \array_module
- *
+ *
* \brief Template functor to compute the cosine of a scalar
*
* \sa class CwiseUnaryOp, MatrixBase::cwiseCos()
@@ -88,7 +116,7 @@ struct ei_functor_traits<ei_scalar_cos_op<Scalar> >
/** \internal
*
* \array_module
- *
+ *
* \brief Template functor to compute the sine of a scalar
*
* \sa class CwiseUnaryOp, MatrixBase::cwiseSin()
@@ -103,7 +131,7 @@ struct ei_functor_traits<ei_scalar_sin_op<Scalar> >
/** \internal
*
* \array_module
- *
+ *
* \brief Template functor to raise a scalar to a power
*
* \sa class CwiseUnaryOp, MatrixBase::cwisePow
@@ -121,7 +149,7 @@ struct ei_functor_traits<ei_scalar_pow_op<Scalar> >
/** \internal
*
* \array_module
- *
+ *
* \brief Template functor to compute the reciprocal of a scalar
*
* \sa class CwiseUnaryOp, MatrixBase::cwiseInverse
diff --git a/Eigen/src/Core/MatrixBase.h b/Eigen/src/Core/MatrixBase.h
index 97c8a39d1..ad0f0b1e8 100644
--- a/Eigen/src/Core/MatrixBase.h
+++ b/Eigen/src/Core/MatrixBase.h
@@ -482,6 +482,9 @@ template<typename Derived> class MatrixBase
/////////// Array module ///////////
+ const Array<Derived> array() const;
+ Array<Derived> array();
+
const CwiseUnaryOp<ei_scalar_sqrt_op<typename ei_traits<Derived>::Scalar>, Derived> cwiseSqrt() const;
const CwiseUnaryOp<ei_scalar_exp_op<typename ei_traits<Derived>::Scalar>, Derived> cwiseExp() const;
const CwiseUnaryOp<ei_scalar_log_op<typename ei_traits<Derived>::Scalar>, Derived> cwiseLog() const;
diff --git a/Eigen/src/Core/util/ForwardDeclarations.h b/Eigen/src/Core/util/ForwardDeclarations.h
index 268b24db0..96949a12a 100644
--- a/Eigen/src/Core/util/ForwardDeclarations.h
+++ b/Eigen/src/Core/util/ForwardDeclarations.h
@@ -55,7 +55,7 @@ template<typename MatrixType> class Map;
template<int Direction, typename UnaryOp, typename MatrixType> class PartialRedux;
template<typename MatrixType, unsigned int Mode> class Part;
template<typename MatrixType, unsigned int Mode> class Extract;
-template<typename Derived, bool HasArrayFlag = int(ei_traits<Derived>::Flags) & ArrayBit> class ArrayBase {};
+template<typename MatrixType> class Array;
template<typename Lhs, typename Rhs> class Cross;
template<typename Scalar> class Quaternion;
template<typename Scalar> class Rotation2D;
diff --git a/Eigen/src/Geometry/Transform.h b/Eigen/src/Geometry/Transform.h
index 0b5b2a3a0..5d092951e 100644
--- a/Eigen/src/Geometry/Transform.h
+++ b/Eigen/src/Geometry/Transform.h
@@ -25,6 +25,16 @@
#ifndef EIGEN_TRANSFORM_H
#define EIGEN_TRANSFORM_H
+// Note that we have to pass Dim and HDim because it is not allowed to use a template
+// parameter to define a template specialization. To be more precise, in the following
+// specializations, it is not allowed to use Dim+1 instead of HDim.
+template< typename Other,
+ int Dim,
+ int HDim,
+ int OtherRows=Other::RowsAtCompileTime,
+ int OtherCols=Other::ColsAtCompileTime>
+struct ei_transform_product_impl;
+
/** \class Transform
*
* \brief Represents an homogeneous transformation in a N dimensional space
@@ -57,52 +67,6 @@ protected:
MatrixType m_matrix;
- template<typename Other,
- int OtherRows=Other::RowsAtCompileTime,
- int OtherCols=Other::ColsAtCompileTime>
- struct ei_transform_product_impl;
-
- // FIXME these specializations of ei_transform_product_impl does not work with gcc 3.3 and 3.4 because
- // Dim depends on a template parameter. Replacing Dim by 3 (for the 3D case) works.
-
- // note that these specializations have to be defined here,
- // otherwise some compilers (at least ICC and NVCC) complain about
- // the use of Dim in the specialization parameters.
- template<typename Other>
- struct ei_transform_product_impl<Other,Dim+1,Dim+1>
- {
- typedef typename Transform<Scalar,Dim>::MatrixType MatrixType;
- typedef typename ProductReturnType<MatrixType,Other>::Type ResultType;
- static ResultType run(const Transform<Scalar,Dim>& tr, const Other& other)
- { return tr.matrix() * other; }
- };
-
- template<typename Other>
- struct ei_transform_product_impl<Other,Dim+1,1>
- {
- typedef typename Transform<Scalar,Dim>::MatrixType MatrixType;
- typedef typename ProductReturnType<MatrixType,Other>::Type ResultType;
- static ResultType run(const Transform<Scalar,Dim>& tr, const Other& other)
- { return tr.matrix() * other; }
- };
-
- template<typename Other>
- struct ei_transform_product_impl<Other,Dim,1>
- {
- typedef typename Transform<Scalar,Dim>::AffineMatrixRef MatrixType;
- typedef const CwiseUnaryOp<
- ei_scalar_multiple_op<Scalar>,
- NestByValue<CwiseBinaryOp<
- ei_scalar_sum_op<Scalar>,
- NestByValue<typename ProductReturnType<NestByValue<MatrixType>,Other>::Type >,
- NestByValue<typename Transform<Scalar,Dim>::VectorRef> > >
- > ResultType;
- // FIXME shall we offer an optimized version when the last row is know to be 0,0...,0,1 ?
- static ResultType run(const Transform<Scalar,Dim>& tr, const Other& other)
- { return ((tr.affine().nestByValue() * other).nestByValue() + tr.translation().nestByValue()).nestByValue()
- * (Scalar(1) / ( (tr.matrix().template block<1,Dim>(Dim,0) * other).coeff(0) + tr.matrix().coeff(Dim,Dim))); }
- };
-
public:
/** Default constructor without initialization of the coefficients. */
@@ -144,7 +108,7 @@ public:
inline VectorRef translation() { return m_matrix.template block<Dim,1>(0,Dim); }
template<typename OtherDerived>
- const typename ei_transform_product_impl<OtherDerived>::ResultType
+ const typename ei_transform_product_impl<OtherDerived,_Dim,_Dim+1>::ResultType
operator * (const MatrixBase<OtherDerived> &other) const;
/** Contatenates two transformations */
@@ -225,12 +189,17 @@ QMatrix Transform<Scalar,Dim>::toQMatrix(void) const
}
#endif
+/** \returns an expression of the product between the transform \c *this and a matrix expression \a other
+ *
+ * The right hand side \a other might be a vector of size Dim, an homogeneous vector of size Dim+1
+ * or a transformation matrix of size Dim+1 x Dim+1.
+ */
template<typename Scalar, int Dim>
template<typename OtherDerived>
-const typename Transform<Scalar,Dim>::template ei_transform_product_impl<OtherDerived>::ResultType
+const typename ei_transform_product_impl<OtherDerived,Dim,Dim+1>::ResultType
Transform<Scalar,Dim>::operator*(const MatrixBase<OtherDerived> &other) const
{
- return ei_transform_product_impl<OtherDerived>::run(*this,other.derived());
+ return ei_transform_product_impl<OtherDerived,Dim,HDim>::run(*this,other.derived());
}
/** Applies on the right the non uniform scale transformation represented
@@ -408,4 +377,47 @@ Transform<Scalar,Dim>::fromPositionOrientationScale(const MatrixBase<PositionDer
return *this;
}
+/***********************************
+*** Specializations of operator* ***
+***********************************/
+
+template<typename Other, int Dim, int HDim>
+struct ei_transform_product_impl<Other,Dim,HDim, HDim,HDim>
+{
+ typedef Transform<typename Other::Scalar,Dim> Transform;
+ typedef typename Transform::MatrixType MatrixType;
+ typedef typename ProductReturnType<MatrixType,Other>::Type ResultType;
+ static ResultType run(const Transform& tr, const Other& other)
+ { return tr.matrix() * other; }
+};
+
+template<typename Other, int Dim, int HDim>
+struct ei_transform_product_impl<Other,Dim,HDim, HDim,1>
+{
+ typedef Transform<typename Other::Scalar,Dim> Transform;
+ typedef typename Transform::MatrixType MatrixType;
+ typedef typename ProductReturnType<MatrixType,Other>::Type ResultType;
+ static ResultType run(const Transform& tr, const Other& other)
+ { return tr.matrix() * other; }
+};
+
+template<typename Other, int Dim, int HDim>
+struct ei_transform_product_impl<Other,Dim,HDim, Dim,1>
+{
+ typedef typename Other::Scalar Scalar;
+ typedef Transform<Scalar,Dim> Transform;
+ typedef typename Transform::AffineMatrixRef MatrixType;
+ typedef const CwiseUnaryOp<
+ ei_scalar_multiple_op<Scalar>,
+ NestByValue<CwiseBinaryOp<
+ ei_scalar_sum_op<Scalar>,
+ NestByValue<typename ProductReturnType<NestByValue<MatrixType>,Other>::Type >,
+ NestByValue<typename Transform::VectorRef> > >
+ > ResultType;
+ // FIXME shall we offer an optimized version when the last row is known to be 0,0...,0,1 ?
+ static ResultType run(const Transform& tr, const Other& other)
+ { return ((tr.affine().nestByValue() * other).nestByValue() + tr.translation().nestByValue()).nestByValue()
+ * (Scalar(1) / ( (tr.matrix().template block<1,Dim>(Dim,0) * other).coeff(0) + tr.matrix().coeff(Dim,Dim))); }
+};
+
#endif // EIGEN_TRANSFORM_H
diff --git a/test/CMakeLists.txt b/test/CMakeLists.txt
index 8dbd3a6af..f89a07e3b 100644
--- a/test/CMakeLists.txt
+++ b/test/CMakeLists.txt
@@ -86,6 +86,7 @@ EI_ADD_TEST(submatrices)
EI_ADD_TEST(miscmatrices)
EI_ADD_TEST(smallvectors)
EI_ADD_TEST(map)
+EI_ADD_TEST(array)
EI_ADD_TEST(triangular)
EI_ADD_TEST(cholesky)
# EI_ADD_TEST(determinant)
diff --git a/test/linearstructure.cpp b/test/linearstructure.cpp
index 02b77c0f0..4c2bdf620 100644
--- a/test/linearstructure.cpp
+++ b/test/linearstructure.cpp
@@ -94,6 +94,7 @@ void test_linearstructure()
{
for(int i = 0; i < g_repeat; i++) {
CALL_SUBTEST( linearStructure(Matrix<float, 1, 1>()) );
+ CALL_SUBTEST( linearStructure(Matrix2f()) );
CALL_SUBTEST( linearStructure(Matrix4d()) );
CALL_SUBTEST( linearStructure(MatrixXcf(3, 3)) );
CALL_SUBTEST( linearStructure(MatrixXf(8, 12)) );