diff options
author | Gael Guennebaud <g.gael@free.fr> | 2008-07-20 15:18:54 +0000 |
---|---|---|
committer | Gael Guennebaud <g.gael@free.fr> | 2008-07-20 15:18:54 +0000 |
commit | ce425d92f1c1d199c8428c46b4c4ea3be81b137b (patch) | |
tree | 582521e6e2f5bd871902e0a3f297acae38695936 | |
parent | 269f683902231020a910fd4e2c1b74554183e2c8 (diff) |
Various documentation improvements, in particualr in Cholesky and Geometry module.
Added doxygen groups for Matrix typedefs and the Geometry module
-rw-r--r-- | Eigen/Geometry | 11 | ||||
-rw-r--r-- | Eigen/src/Cholesky/Cholesky.h | 12 | ||||
-rw-r--r-- | Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h | 18 | ||||
-rwxr-xr-x | Eigen/src/Core/InverseProduct.h | 4 | ||||
-rw-r--r-- | Eigen/src/Core/Matrix.h | 17 | ||||
-rw-r--r-- | Eigen/src/Core/Part.h | 8 | ||||
-rw-r--r-- | Eigen/src/Core/Visitor.h | 3 | ||||
-rw-r--r-- | Eigen/src/Core/util/Constants.h | 4 | ||||
-rw-r--r-- | Eigen/src/Geometry/AngleAxis.h | 11 | ||||
-rw-r--r-- | Eigen/src/Geometry/Quaternion.h | 59 | ||||
-rw-r--r-- | Eigen/src/Geometry/Rotation.h | 20 | ||||
-rw-r--r-- | Eigen/src/Geometry/Transform.h | 95 | ||||
-rw-r--r-- | doc/Doxyfile.in | 4 | ||||
-rw-r--r-- | doc/snippets/AngleAxis_mimic_euler.cpp | 2 |
14 files changed, 179 insertions, 89 deletions
diff --git a/Eigen/Geometry b/Eigen/Geometry index 429dc2ac5..5e73ce992 100644 --- a/Eigen/Geometry +++ b/Eigen/Geometry @@ -5,7 +5,16 @@ namespace Eigen { -/** \defgroup Geometry */ +/** \defgroup Geometry + * This module provides support for: + * - fixed-size homogeneous transformations + * - 2D and 3D rotations + * - \ref MatrixBase::cross() "cross product" + * + * \code + * #include <Eigen/Geometry> + * \endcode + */ #include "src/Geometry/Cross.h" #include "src/Geometry/Quaternion.h" diff --git a/Eigen/src/Cholesky/Cholesky.h b/Eigen/src/Cholesky/Cholesky.h index c1f05d768..dd4fc6e38 100644 --- a/Eigen/src/Cholesky/Cholesky.h +++ b/Eigen/src/Cholesky/Cholesky.h @@ -42,7 +42,7 @@ * Note that during the decomposition, only the upper triangular part of A is considered. Therefore, * the strict lower part does not have to store correct values. * - * \sa class CholeskyWithoutSquareRoot + * \sa MatrixBase::cholesky(), class CholeskyWithoutSquareRoot */ template<typename MatrixType> class Cholesky { @@ -107,20 +107,22 @@ void Cholesky<MatrixType>::compute(const MatrixType& a) } } -/** \returns the solution of A x = \a b using the current decomposition of A. - * In other words, it returns \code A^-1 b \endcode computing - * \code L^-* L^1 b \endcode from right to left. +/** \returns the solution of \f$ A x = b \f$ using the current decomposition of A. + * In other words, it returns \f$ A^{-1} b \f$ computing + * \f$ {L^{*}}^{-1} L^{-1} b \f$ from right to left. + * \param b the column vector \f$ b \f$, which can also be a matrix. * * Example: \include Cholesky_solve.cpp * Output: \verbinclude Cholesky_solve.out * + * \sa MatrixBase::cholesky(), CholeskyWithoutSquareRoot::solve() */ template<typename MatrixType> template<typename Derived> typename Derived::Eval Cholesky<MatrixType>::solve(const MatrixBase<Derived> &b) const { const int size = m_matrix.rows(); - ei_assert(size==b.size()); + ei_assert(size==b.rows()); return m_matrix.adjoint().template extract<Upper>().inverseProduct(matrixL().inverseProduct(b)); } diff --git a/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h b/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h index 2d85f78db..2572b88a2 100644 --- a/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h +++ b/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h @@ -32,8 +32,8 @@ * \param MatrixType the type of the matrix of which we are computing the Cholesky decomposition * * This class performs a Cholesky decomposition without square root of a symmetric, positive definite - * matrix A such that A = L D L^* = U^* D U, where L is lower triangular with a unit diagonal and D is a diagonal - * matrix. + * matrix A such that A = L D L^* = U^* D U, where L is lower triangular with a unit diagonal + * and D is a diagonal matrix. * * Compared to a standard Cholesky decomposition, avoiding the square roots allows for faster and more * stable computation. @@ -41,7 +41,7 @@ * Note that during the decomposition, only the upper triangular part of A is considered. Therefore, * the strict lower part does not have to store correct values. * - * \sa class Cholesky + * \sa MatrixBase::choleskyNoSqrt(), class Cholesky */ template<typename MatrixType> class CholeskyWithoutSquareRoot { @@ -123,19 +123,23 @@ void CholeskyWithoutSquareRoot<MatrixType>::compute(const MatrixType& a) /** \returns the solution of \f$ A x = b \f$ using the current decomposition of A. * In other words, it returns \f$ A^{-1} b \f$ computing * \f$ {L^{*}}^{-1} D^{-1} L^{-1} b \f$ from right to left. - * \param vecB the vector \f$ b \f$ (or an array of vectors) + * \param b the column vector \f$ b \f$, which can also be a matrix. + * + * See Cholesky::solve() for a example. + * + * \sa MatrixBase::choleskyNoSqrt() */ template<typename MatrixType> template<typename Derived> -typename Derived::Eval CholeskyWithoutSquareRoot<MatrixType>::solve(const MatrixBase<Derived> &vecB) const +typename Derived::Eval CholeskyWithoutSquareRoot<MatrixType>::solve(const MatrixBase<Derived> &b) const { const int size = m_matrix.rows(); - ei_assert(size==vecB.size()); + ei_assert(size==b.rows()); return m_matrix.adjoint().template extract<UnitUpper>() .inverseProduct( (matrixL() - .inverseProduct(vecB)) + .inverseProduct(b)) .cwise()/m_matrix.diagonal() ); } diff --git a/Eigen/src/Core/InverseProduct.h b/Eigen/src/Core/InverseProduct.h index 57dbf7509..0ee54a3fb 100755 --- a/Eigen/src/Core/InverseProduct.h +++ b/Eigen/src/Core/InverseProduct.h @@ -72,9 +72,9 @@ void MatrixBase<Derived>::inverseProductInPlace(MatrixBase<OtherDerived>& other) } } -/** \returns the product of the inverse of \c *this with \a other. +/** \returns the product of the inverse of \c *this with \a other, \a *this being triangular. * - * This function computes the inverse-matrix matrix product inverse(\c*this) * \a other + * This function computes the inverse-matrix matrix product inverse(\c *this) * \a other * It works as a forward (resp. backward) substitution if \c *this is an upper (resp. lower) * triangular matrix. * diff --git a/Eigen/src/Core/Matrix.h b/Eigen/src/Core/Matrix.h index 92988b725..85a984920 100644 --- a/Eigen/src/Core/Matrix.h +++ b/Eigen/src/Core/Matrix.h @@ -71,6 +71,8 @@ * \li \c VectorXf is a typedef for \c Matrix<float,Dynamic,1> * \li \c RowVector3i is a typedef for \c Matrix<int,1,3> * + * See \ref matrixtypedefs for an explicit list of all matrix typedefs. + * * Of course these typedefs do not exhaust all the possibilities offered by the Matrix class * template, they only address some of the most common cases. For instance, if you want a * fixed-size matrix with 3 rows and 5 columns, there is no typedef for that, so you should use @@ -355,9 +357,18 @@ class Matrix : public MatrixBase<Matrix<_Scalar, _Rows, _Cols, _MaxRows, _MaxCol } }; -#define EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, Size, SizeSuffix) \ -typedef Matrix<Type, Size, Size> Matrix##SizeSuffix##TypeSuffix; \ -typedef Matrix<Type, Size, 1> Vector##SizeSuffix##TypeSuffix; \ +/** \defgroup matrixtypedefs Global matrix typedefs + * Eigen defines several typedef shortcuts for most common matrix types. + * Here is the explicit list. + * \sa class Matrix + */ + +#define EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, Size, SizeSuffix) \ +/** \ingroup matrixtypedefs */ \ +typedef Matrix<Type, Size, Size> Matrix##SizeSuffix##TypeSuffix; \ +/** \ingroup matrixtypedefs */ \ +typedef Matrix<Type, Size, 1> Vector##SizeSuffix##TypeSuffix; \ +/** \ingroup matrixtypedefs */ \ typedef Matrix<Type, 1, Size> RowVector##SizeSuffix##TypeSuffix; #define EIGEN_MAKE_TYPEDEFS_ALL_SIZES(Type, TypeSuffix) \ diff --git a/Eigen/src/Core/Part.h b/Eigen/src/Core/Part.h index eb8dcbba7..a0553923a 100644 --- a/Eigen/src/Core/Part.h +++ b/Eigen/src/Core/Part.h @@ -106,9 +106,9 @@ struct ei_part_assignment_impl if(Mode == SelfAdjoint) { if(row == col) - dst.coeffRef(row, col) = ei_real(src.coeff(row, col)); + dst.coeffRef(row, col) = ei_real(src.coeff(row, col)); else if(row < col) - dst.coeffRef(col, row) = ei_conj(dst.coeffRef(row, col) = src.coeff(row, col)); + dst.coeffRef(col, row) = ei_conj(dst.coeffRef(row, col) = src.coeff(row, col)); } else { @@ -116,7 +116,7 @@ struct ei_part_assignment_impl || (Mode == Lower && row >= col) || (Mode == StrictlyUpper && row < col) || (Mode == StrictlyLower && row > col)) - dst.coeffRef(row, col) = src.coeff(row, col); + dst.coeffRef(row, col) = src.coeff(row, col); } } }; @@ -262,6 +262,8 @@ inline void Part<MatrixType, Mode>::setRandom() * The \a Mode parameter can have the following values: \c Upper, \c StrictlyUpper, \c Lower, * \c StrictlyLower, \c SelfAdjoint. * + * \addexample PartExample \label How to write to a triangular part of a matrix + * * Example: \include MatrixBase_part.cpp * Output: \verbinclude MatrixBase_part.out * diff --git a/Eigen/src/Core/Visitor.h b/Eigen/src/Core/Visitor.h index 18e90ca62..041aa9445 100644 --- a/Eigen/src/Core/Visitor.h +++ b/Eigen/src/Core/Visitor.h @@ -76,6 +76,9 @@ struct ei_visitor_impl<Visitor, Derived, Dynamic> * }; * \endcode * + * \note compared to one or two \em for \em loops, visitors offer automatic + * unrolling for small fixed size matrix. + * * \sa minCoeff(int*,int*), maxCoeff(int*,int*), MatrixBase::redux() */ template<typename Derived> diff --git a/Eigen/src/Core/util/Constants.h b/Eigen/src/Core/util/Constants.h index 7e2d37dff..24c653e2e 100644 --- a/Eigen/src/Core/util/Constants.h +++ b/Eigen/src/Core/util/Constants.h @@ -28,9 +28,7 @@ const int Dynamic = 10000; -/** \defgroup flags */ - -/** \name flags +/** \defgroup flags * * These are the possible bits which can be OR'ed to constitute the flags of a matrix or * expression. diff --git a/Eigen/src/Geometry/AngleAxis.h b/Eigen/src/Geometry/AngleAxis.h index 647e07513..ca53140fb 100644 --- a/Eigen/src/Geometry/AngleAxis.h +++ b/Eigen/src/Geometry/AngleAxis.h @@ -64,10 +64,15 @@ protected: public: + /** Default constructor without initialization. */ AngleAxis() {} + /** Constructs and initialize the angle-axis rotation from an \a angle in radian + * and an \a axis which must be normalized. */ template<typename Derived> inline AngleAxis(Scalar angle, const MatrixBase<Derived>& axis) : m_axis(axis), m_angle(angle) {} + /** Constructs and initialize the angle-axis rotation from a quaternion \a q. */ inline AngleAxis(const QuaternionType& q) { *this = q; } + /** Constructs and initialize the angle-axis rotation from a 3x3 rotation matrix. */ template<typename Derived> inline AngleAxis(const MatrixBase<Derived>& m) { *this = m; } @@ -77,6 +82,8 @@ public: const Vector3& axis() const { return m_axis; } Vector3& axis() { return m_axis; } + /** Automatic conversion to a 3x3 rotation matrix. + * \sa toRotationMatrix() */ operator Matrix3 () const { return toRotationMatrix(); } inline QuaternionType operator* (const AngleAxis& other) const @@ -105,7 +112,11 @@ public: Matrix3 toRotationMatrix(void) const; }; +/** \ingroup Geometry + * single precision angle-axis type */ typedef AngleAxis<float> AngleAxisf; +/** \ingroup Geometry + * double precision angle-axis type */ typedef AngleAxis<double> AngleAxisd; /** Set \c *this from a quaternion. diff --git a/Eigen/src/Geometry/Quaternion.h b/Eigen/src/Geometry/Quaternion.h index 1d31e2f96..ab6994b4b 100644 --- a/Eigen/src/Geometry/Quaternion.h +++ b/Eigen/src/Geometry/Quaternion.h @@ -38,13 +38,12 @@ struct ei_quaternion_assign_impl; * * \param _Scalar the scalar type, i.e., the type of the coefficients * - * This class represents a quaternion that is a convenient representation of - * orientations and rotations of objects in three dimensions. Compared to other - * representations like Euler angles or 3x3 matrices, quatertions offer the - * following advantages: - * \li \c compact storage (4 scalars) - * \li \c efficient to compose (28 flops), - * \li \c stable spherical interpolation + * This class represents a quaternion \f$ w+xi+yj+zk \f$ that is a convenient representation of + * orientations and rotations of objects in three dimensions. Compared to other representations + * like Euler angles or 3x3 matrices, quatertions offer the following advantages: + * \li \b compact storage (4 scalars) + * \li \b efficient to compose (28 flops), + * \li \b stable spherical interpolation * * The following two typedefs are provided for convenience: * \li \c Quaternionf for \c float @@ -63,18 +62,29 @@ public: /** the scalar type of the coefficients */ typedef _Scalar Scalar; + /** the type of a 3D vector */ typedef Matrix<Scalar,3,1> Vector3; + /** the equivalent rotation matrix type */ typedef Matrix<Scalar,3,3> Matrix3; + /** the equivalent angle-axis type */ typedef AngleAxis<Scalar> AngleAxisType; + /** \returns the \c x coefficient */ inline Scalar x() const { return m_coeffs.coeff(0); } + /** \returns the \c y coefficient */ inline Scalar y() const { return m_coeffs.coeff(1); } + /** \returns the \c z coefficient */ inline Scalar z() const { return m_coeffs.coeff(2); } + /** \returns the \c w coefficient */ inline Scalar w() const { return m_coeffs.coeff(3); } + /** \returns a reference to the \c x coefficient */ inline Scalar& x() { return m_coeffs.coeffRef(0); } + /** \returns a reference to the \c y coefficient */ inline Scalar& y() { return m_coeffs.coeffRef(1); } + /** \returns a reference to the \c z coefficient */ inline Scalar& z() { return m_coeffs.coeffRef(2); } + /** \returns a reference to the \c w coefficient */ inline Scalar& w() { return m_coeffs.coeffRef(3); } /** \returns a read-only vector expression of the imaginary part (x,y,z) */ @@ -83,25 +93,33 @@ public: /** \returns a vector expression of the imaginary part (x,y,z) */ inline Block<Coefficients,3,1> vec() { return m_coeffs.template start<3>(); } - /** \returns a read-only vector expression of the coefficients */ + /** \returns a read-only vector expression of the coefficients (x,y,z,w) */ inline const Coefficients& coeffs() const { return m_coeffs; } - /** \returns a vector expression of the coefficients */ + /** \returns a vector expression of the coefficients (x,y,z,w) */ inline Coefficients& coeffs() { return m_coeffs; } + /** Default constructor and initializing an identity quaternion. */ + inline Quaternion() + { m_coeffs << 0, 0, 0, 1; } + + /** Constructs and initializes the quaternion \f$ w+xi+yj+zk \f$ from + * its four coefficients \a w, \a x, \a y and \a z. + */ // FIXME what is the prefered order: w x,y,z or x,y,z,w ? - inline Quaternion(Scalar w = 1.0, Scalar x = 0.0, Scalar y = 0.0, Scalar z = 0.0) - { - m_coeffs.coeffRef(0) = x; - m_coeffs.coeffRef(1) = y; - m_coeffs.coeffRef(2) = z; - m_coeffs.coeffRef(3) = w; - } + inline Quaternion(Scalar w, Scalar x, Scalar y, Scalar z) + { m_coeffs << x, y, z, w; } /** Copy constructor */ inline Quaternion(const Quaternion& other) { m_coeffs = other.m_coeffs; } + /** Constructs and initializes a quaternion from the angle-axis \a aa */ explicit inline Quaternion(const AngleAxisType& aa) { *this = aa; } + /** Constructs and initializes a quaternion from either: + * - a rotation matrix expression, + * - a 4D vector expression representing quaternion coefficients. + * \sa operator=(MatrixBase<Derived>) + */ template<typename Derived> explicit inline Quaternion(const MatrixBase<Derived>& other) { *this = other; } @@ -110,6 +128,7 @@ public: template<typename Derived> Quaternion& operator=(const MatrixBase<Derived>& m); + /** Automatic conversion to a rotation matrix. */ operator Matrix3 () const { return toRotationMatrix(); } /** \returns a quaternion representing an identity rotation @@ -149,7 +168,11 @@ public: }; +/** \ingroup Geometry + * single precision quaternion type */ typedef Quaternion<float> Quaternionf; +/** \ingroup Geometry + * double precision quaternion type */ typedef Quaternion<double> Quaterniond; /** \returns the concatenation of two rotations as a quaternion-quaternion product */ @@ -165,6 +188,7 @@ inline Quaternion<Scalar> Quaternion<Scalar>::operator* (const Quaternion& other ); } +/** \sa operator*(Quaternion) */ template <typename Scalar> inline Quaternion<Scalar>& Quaternion<Scalar>::operator*= (const Quaternion& other) { @@ -200,8 +224,7 @@ inline Quaternion<Scalar>& Quaternion<Scalar>::operator=(const Quaternion& other return *this; } -/** Set \c *this from an angle-axis \a aa - * and returns a reference to \c *this +/** Set \c *this from an angle-axis \a aa and returns a reference to \c *this */ template<typename Scalar> inline Quaternion<Scalar>& Quaternion<Scalar>::operator=(const AngleAxisType& aa) diff --git a/Eigen/src/Geometry/Rotation.h b/Eigen/src/Geometry/Rotation.h index e0855b739..a98ed061a 100644 --- a/Eigen/src/Geometry/Rotation.h +++ b/Eigen/src/Geometry/Rotation.h @@ -28,7 +28,7 @@ // this file aims to contains the various representations of rotation/orientation // in 2D and 3D space excepted Matrix and Quaternion. -/** \geometry_module +/** \internal * * \class ToRotationMatrix * @@ -103,7 +103,7 @@ struct ToRotationMatrix<Scalar, Dim, MatrixBase<OtherDerived> > } }; -/** \geometry_module +/** \geometry_module \ingroup Geometry * * \class Rotation2D * @@ -111,10 +111,10 @@ struct ToRotationMatrix<Scalar, Dim, MatrixBase<OtherDerived> > * * \param _Scalar the scalar type, i.e., the type of the coefficients * - * This class is equivalent to a single scalar representing the rotation angle - * in radian with some additional features such as the conversion from/to - * rotation matrix. Moreover this class aims to provide a similar interface - * to Quaternion in order to facilitate the writing of generic algorithm + * This class is equivalent to a single scalar representing a counter clock wise rotation + * as a single angle in radian. It provides some additional features such as the automatic + * conversion from/to a 2x2 rotation matrix. Moreover this class aims to provide a similar + * interface to Quaternion in order to facilitate the writing of generic algorithm * dealing with rotations. * * \sa class Quaternion, class Transform @@ -134,16 +134,22 @@ protected: public: + /** Construct a 2D counter clock wise rotation from the angle \a a in radian. */ inline Rotation2D(Scalar a) : m_angle(a) {} inline operator Scalar& () { return m_angle; } inline operator Scalar () const { return m_angle; } + /** Automatic convertion to a 2D rotation matrix. + * \sa toRotationMatrix() + */ + inline operator Matrix2() const { return toRotationMatrix(); } + template<typename Derived> Rotation2D& fromRotationMatrix(const MatrixBase<Derived>& m); Matrix2 toRotationMatrix(void) const; /** \returns the spherical interpolation between \c *this and \a other using - * parameter \a t. It is equivalent to a linear interpolation. + * parameter \a t. It is in fact equivalent to a linear interpolation. */ inline Rotation2D slerp(Scalar t, const Rotation2D& other) const { return m_angle * (1-t) + t * other; } diff --git a/Eigen/src/Geometry/Transform.h b/Eigen/src/Geometry/Transform.h index 7669838de..73e268e13 100644 --- a/Eigen/src/Geometry/Transform.h +++ b/Eigen/src/Geometry/Transform.h @@ -50,20 +50,29 @@ struct ei_transform_product_impl; * Conversion methods from/to Qt's QMatrix are available if the preprocessor token * EIGEN_QT_SUPPORT is defined. * + * \sa class Matrix, class Quaternion */ template<typename _Scalar, int _Dim> class Transform { public: - enum { Dim = _Dim, HDim = _Dim+1 }; + enum { + Dim = _Dim, ///< space dimension in which the transformation holds + HDim = _Dim+1 ///< size of a respective homogeneous vector + }; /** the scalar type of the coefficients */ typedef _Scalar Scalar; + /** type of the matrix used to represent the transformation */ typedef Matrix<Scalar,HDim,HDim> MatrixType; + /** type of the matrix used to represent the affine part of the transformation */ typedef Matrix<Scalar,Dim,Dim> AffineMatrixType; - typedef Block<MatrixType,Dim,Dim> AffineMatrixRef; + /** type of read/write reference to the affine part of the transformation */ + typedef Block<MatrixType,Dim,Dim> AffinePart; + /** type of a vector */ typedef Matrix<Scalar,Dim,1> VectorType; - typedef Block<MatrixType,Dim,1> VectorRef; + /** type of a read/write reference to the translation part of the rotation */ + typedef Block<MatrixType,Dim,1> TranslationPart; protected: @@ -80,10 +89,12 @@ public: inline Transform& operator=(const Transform& other) { m_matrix = other.m_matrix; return *this; } + /** Constructs and initializes a transformation from a (Dim+1)^2 matrix. */ template<typename OtherDerived> inline explicit Transform(const MatrixBase<OtherDerived>& other) { m_matrix = other; } + /** Set \c *this from a (Dim+1)^2 matrix. */ template<typename OtherDerived> inline Transform& operator=(const MatrixBase<OtherDerived>& other) { m_matrix = other; return *this; } @@ -100,14 +111,14 @@ public: inline MatrixType& matrix() { return m_matrix; } /** \returns a read-only expression of the affine (linear) part of the transformation */ - inline const AffineMatrixRef affine() const { return m_matrix.template block<Dim,Dim>(0,0); } + inline const AffinePart affine() const { return m_matrix.template block<Dim,Dim>(0,0); } /** \returns a writable expression of the affine (linear) part of the transformation */ - inline AffineMatrixRef affine() { return m_matrix.template block<Dim,Dim>(0,0); } + inline AffinePart affine() { return m_matrix.template block<Dim,Dim>(0,0); } /** \returns a read-only expression of the translation vector of the transformation */ - inline const VectorRef translation() const { return m_matrix.template block<Dim,1>(0,Dim); } + inline const TranslationPart translation() const { return m_matrix.template block<Dim,1>(0,Dim); } /** \returns a writable expression of the translation vector of the transformation */ - inline VectorRef translation() { return m_matrix.template block<Dim,1>(0,Dim); } + inline TranslationPart translation() { return m_matrix.template block<Dim,1>(0,Dim); } template<typename OtherDerived> const typename ei_transform_product_impl<OtherDerived,_Dim,_Dim+1>::ResultType @@ -118,6 +129,7 @@ public: operator * (const Transform& other) const { return m_matrix * other.matrix(); } + /** \sa MatrixBase::setIdentity() */ void setIdentity() { m_matrix.setIdentity(); } template<typename OtherDerived> @@ -138,10 +150,7 @@ public: template<typename RotationType> Transform& prerotate(const RotationType& rotation); - template<typename OtherDerived> Transform& shear(Scalar sx, Scalar sy); - - template<typename OtherDerived> Transform& preshear(Scalar sx, Scalar sy); AffineMatrixType extractRotation() const; @@ -151,6 +160,7 @@ public: Transform& fromPositionOrientationScale(const MatrixBase<PositionDerived> &position, const OrientationType& orientation, const MatrixBase<ScaleDerived> &scale); + /** \sa MatrixBase::inverse() */ const Inverse<MatrixType, false> inverse() const { return m_matrix.inverse(); } @@ -158,8 +168,19 @@ protected: }; +/** \ingroup Geometry */ +typedef Transform<float,2> Transform2f; +/** \ingroup Geometry */ +typedef Transform<float,3> Transform3f; +/** \ingroup Geometry */ +typedef Transform<double,2> Transform2d; +/** \ingroup Geometry */ +typedef Transform<double,3> Transform3d; + #ifdef EIGEN_QT_SUPPORT /** Initialises \c *this from a QMatrix assuming the dimension is 2. + * + * This function is available only if the token EIGEN_QT_SUPPORT is defined. */ template<typename Scalar, int Dim> Transform<Scalar,Dim>::Transform(const QMatrix& other) @@ -168,6 +189,8 @@ Transform<Scalar,Dim>::Transform(const QMatrix& other) } /** Set \c *this from a QMatrix assuming the dimension is 2. + * + * This function is available only if the token EIGEN_QT_SUPPORT is defined. */ template<typename Scalar, int Dim> Transform<Scalar,Dim>& Transform<Scalar,Dim>::operator=(const QMatrix& other) @@ -180,21 +203,25 @@ Transform<Scalar,Dim>& Transform<Scalar,Dim>::operator=(const QMatrix& other) } /** \returns a QMatrix from \c *this assuming the dimension is 2. + * + * This function is available only if the token EIGEN_QT_SUPPORT is defined. */ template<typename Scalar, int Dim> QMatrix Transform<Scalar,Dim>::toQMatrix(void) const { EIGEN_STATIC_ASSERT(Dim==2, you_did_a_programming_error); - return QMatrix( other.coeffRef(0,0), other.coeffRef(1,0), - other.coeffRef(0,1), other.coeffRef(1,1), - other.coeffRef(0,2), other.coeffRef(1,2), + return QMatrix(other.coeffRef(0,0), other.coeffRef(1,0), + other.coeffRef(0,1), other.coeffRef(1,1), + other.coeffRef(0,2), other.coeffRef(1,2)); } #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. + * The right hand side \a other might be either: + * \li a vector of size Dim, + * \li an homogeneous vector of size Dim+1, + * \li a transformation matrix of size Dim+1 x Dim+1. */ template<typename Scalar, int Dim> template<typename OtherDerived> @@ -213,8 +240,7 @@ template<typename OtherDerived> Transform<Scalar,Dim>& Transform<Scalar,Dim>::scale(const MatrixBase<OtherDerived> &other) { - EIGEN_STATIC_ASSERT(int(OtherDerived::IsVectorAtCompileTime) - && int(OtherDerived::SizeAtCompileTime)==int(Dim), you_did_a_programming_error); + EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,Dim); affine() = (affine() * other.asDiagonal()).lazy(); return *this; } @@ -228,8 +254,7 @@ template<typename OtherDerived> Transform<Scalar,Dim>& Transform<Scalar,Dim>::prescale(const MatrixBase<OtherDerived> &other) { - EIGEN_STATIC_ASSERT(int(OtherDerived::IsVectorAtCompileTime) - && int(OtherDerived::SizeAtCompileTime)==int(Dim), you_did_a_programming_error); + EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,Dim); m_matrix.template block<Dim,HDim>(0,0) = (other.asDiagonal() * m_matrix.template block<Dim,HDim>(0,0)).lazy(); return *this; } @@ -243,8 +268,7 @@ template<typename OtherDerived> Transform<Scalar,Dim>& Transform<Scalar,Dim>::translate(const MatrixBase<OtherDerived> &other) { - EIGEN_STATIC_ASSERT(int(OtherDerived::IsVectorAtCompileTime) - && int(OtherDerived::SizeAtCompileTime)==int(Dim), you_did_a_programming_error); + EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,Dim); translation() += affine() * other; return *this; } @@ -258,8 +282,7 @@ template<typename OtherDerived> Transform<Scalar,Dim>& Transform<Scalar,Dim>::pretranslate(const MatrixBase<OtherDerived> &other) { - EIGEN_STATIC_ASSERT(int(OtherDerived::IsVectorAtCompileTime) - && int(OtherDerived::SizeAtCompileTime)==int(Dim), you_did_a_programming_error); + EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,Dim); translation() += other; return *this; } @@ -273,8 +296,8 @@ Transform<Scalar,Dim>::pretranslate(const MatrixBase<OtherDerived> &other) * Natively supported types includes: * - any scalar (2D), * - a Dim x Dim matrix expression, - * - Quaternion (3D), - * - AngleAxis (3D) + * - a Quaternion (3D), + * - a AngleAxis (3D) * * This mechanism is easily extendable to support user types such as Euler angles, * or a pair of Quaternion for 4D rotations. @@ -293,9 +316,9 @@ Transform<Scalar,Dim>::rotate(const RotationType& rotation) /** Applies on the left the rotation represented by the rotation \a rotation * to \c *this and returns a reference to \c *this. * - * See rotate(RotationType) for further details. + * See rotate() for further details. * - * \sa rotate(RotationType), rotate(Scalar) + * \sa rotate() */ template<typename Scalar, int Dim> template<typename RotationType> @@ -303,7 +326,7 @@ Transform<Scalar,Dim>& Transform<Scalar,Dim>::prerotate(const RotationType& rotation) { m_matrix.template block<Dim,HDim>(0,0) = ToRotationMatrix<Scalar,Dim,RotationType>::convert(rotation) - * m_matrix.template block<Dim,HDim>(0,0); + * m_matrix.template block<Dim,HDim>(0,0); return *this; } @@ -313,12 +336,10 @@ Transform<Scalar,Dim>::prerotate(const RotationType& rotation) * \sa preshear() */ template<typename Scalar, int Dim> -template<typename OtherDerived> Transform<Scalar,Dim>& Transform<Scalar,Dim>::shear(Scalar sx, Scalar sy) { - EIGEN_STATIC_ASSERT(int(OtherDerived::IsVectorAtCompileTime) - && int(OtherDerived::SizeAtCompileTime)==int(Dim) && int(Dim)==2, you_did_a_programming_error); + EIGEN_STATIC_ASSERT(int(Dim)==2, you_did_a_programming_error); VectorType tmp = affine().col(0)*sy + affine().col(1); affine() << affine().col(0) + affine().col(1)*sx, tmp; return *this; @@ -330,18 +351,16 @@ Transform<Scalar,Dim>::shear(Scalar sx, Scalar sy) * \sa shear() */ template<typename Scalar, int Dim> -template<typename OtherDerived> Transform<Scalar,Dim>& Transform<Scalar,Dim>::preshear(Scalar sx, Scalar sy) { - EIGEN_STATIC_ASSERT(int(OtherDerived::IsVectorAtCompileTime) - && int(OtherDerived::SizeAtCompileTime)==int(Dim), you_did_a_programming_error); + EIGEN_STATIC_ASSERT(int(Dim)==2, you_did_a_programming_error); m_matrix.template block<Dim,HDim>(0,0) = AffineMatrixType(1, sx, sy, 1) * m_matrix.template block<Dim,HDim>(0,0); return *this; } /** \returns the rotation part of the transformation using a QR decomposition. - * \sa extractRotationNoShear() + * \sa extractRotationNoShear(), class QR */ template<typename Scalar, int Dim> typename Transform<Scalar,Dim>::AffineMatrixType @@ -408,15 +427,15 @@ struct ei_transform_product_impl<Other,Dim,HDim, Dim,1> { typedef typename Other::Scalar Scalar; typedef Transform<Scalar,Dim> TransformType; - typedef typename TransformType::AffineMatrixRef MatrixType; + typedef typename TransformType::AffinePart 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 TransformType::VectorRef> > > + NestByValue<typename TransformType::TranslationPart> > > > ResultType; - // FIXME shall we offer an optimized version when the last row is known to be 0,0...,0,1 ? + // FIXME should we offer an optimized version when the last row is known to be 0,0...,0,1 ? static ResultType run(const TransformType& 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))); } diff --git a/doc/Doxyfile.in b/doc/Doxyfile.in index 4c5c71fd7..eed3fe600 100644 --- a/doc/Doxyfile.in +++ b/doc/Doxyfile.in @@ -1191,7 +1191,9 @@ PREDEFINED = EIGEN_EMPTY_STRUCT \ # The macro definition that is found in the sources will be used. # Use the PREDEFINED tag if you want to use a different macro definition. -EXPAND_AS_DEFINED = EIGEN_MAKE_SCALAR_OPS +EXPAND_AS_DEFINED = EIGEN_MAKE_SCALAR_OPS \ + EIGEN_MAKE_TYPEDEFS \ + EIGEN_MAKE_TYPEDEFS_ALL_SIZES # If the SKIP_FUNCTION_MACROS tag is set to YES (the default) then # doxygen's preprocessor will remove all function-like macros that are alone diff --git a/doc/snippets/AngleAxis_mimic_euler.cpp b/doc/snippets/AngleAxis_mimic_euler.cpp index be6b8adbe..46ec41aa5 100644 --- a/doc/snippets/AngleAxis_mimic_euler.cpp +++ b/doc/snippets/AngleAxis_mimic_euler.cpp @@ -1,4 +1,4 @@ Matrix3f m = AngleAxisf(0.25*M_PI, Vector3f::UnitX()) * AngleAxisf(0.5*M_PI, Vector3f::UnitY()) * AngleAxisf(0.33*M_PI, Vector3f::UnitZ()); -cout << m << endl; +cout << m << endl << "is unitary: " << m.isUnitary() << endl; |