diff options
author | Gael Guennebaud <g.gael@free.fr> | 2009-02-17 09:53:05 +0000 |
---|---|---|
committer | Gael Guennebaud <g.gael@free.fr> | 2009-02-17 09:53:05 +0000 |
commit | e6f1104b57f19dff773b4f22d26d6aacabd1bdb2 (patch) | |
tree | 5de1c7a6d6eb47d894f119b5bda0107bd531000d /Eigen/src/Geometry | |
parent | 67b4fab4e30a59d9a7e001ef25938d1767371569 (diff) |
* fix Quaternion::setFromTwoVectors (thanks to "benv" from the forum)
* extend PartialRedux::cross() to any matrix sizes with automatic
vectorization when possible
* unit tests: add "geo_" prefix to all unit tests related to the
geometry module and start splitting the big "geometry.cpp" tests to
multiple smaller ones (also include new tests)
Diffstat (limited to 'Eigen/src/Geometry')
-rw-r--r-- | Eigen/src/Geometry/OrthoMethods.h | 38 | ||||
-rw-r--r-- | Eigen/src/Geometry/Quaternion.h | 11 |
2 files changed, 47 insertions, 2 deletions
diff --git a/Eigen/src/Geometry/OrthoMethods.h b/Eigen/src/Geometry/OrthoMethods.h index 047152d0b..e22df114a 100644 --- a/Eigen/src/Geometry/OrthoMethods.h +++ b/Eigen/src/Geometry/OrthoMethods.h @@ -1,7 +1,7 @@ // 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> +// Copyright (C) 2008-2009 Gael Guennebaud <g.gael@free.fr> // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com> // // Eigen is free software; you can redistribute it and/or @@ -50,6 +50,42 @@ MatrixBase<Derived>::cross(const MatrixBase<OtherDerived>& other) const ); } +/** \returns a matrix expression of the cross product of each column or row + * of the referenced expression with the \a other vector. + * + * The referenced matrix must have one dimension equal to 3. + * The result matrix has the same dimensions than the referenced one. + * + * \geometry_module + * + * \sa MatrixBase::cross() */ +template<typename ExpressionType, int Direction> +template<typename OtherDerived> +const typename PartialRedux<ExpressionType,Direction>::CrossReturnType +PartialRedux<ExpressionType,Direction>::cross(const MatrixBase<OtherDerived>& other) const +{ + EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3) + EIGEN_STATIC_ASSERT((ei_is_same_type<Scalar, typename OtherDerived::Scalar>::ret), + YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY) + + CrossReturnType res(_expression().rows(),_expression().cols()); + if(Direction==Vertical) + { + ei_assert(CrossReturnType::RowsAtCompileTime==3 && "the matrix must have exactly 3 rows"); + res.row(0) = _expression().row(1) * other.coeff(2) - _expression().row(2) * other.coeff(1); + res.row(1) = _expression().row(2) * other.coeff(0) - _expression().row(0) * other.coeff(2); + res.row(2) = _expression().row(0) * other.coeff(1) - _expression().row(1) * other.coeff(0); + } + else + { + ei_assert(CrossReturnType::ColsAtCompileTime==3 && "the matrix must have exactly 3 columns"); + res.col(0) = _expression().col(1) * other.coeff(2) - _expression().col(2) * other.coeff(1); + res.col(1) = _expression().col(2) * other.coeff(0) - _expression().col(0) * other.coeff(2); + res.col(2) = _expression().col(0) * other.coeff(1) - _expression().col(1) * other.coeff(0); + } + return res; +} + template<typename Derived, int Size = Derived::SizeAtCompileTime> struct ei_unitOrthogonal_selector { diff --git a/Eigen/src/Geometry/Quaternion.h b/Eigen/src/Geometry/Quaternion.h index e72ce814f..0305fe176 100644 --- a/Eigen/src/Geometry/Quaternion.h +++ b/Eigen/src/Geometry/Quaternion.h @@ -345,7 +345,6 @@ inline Quaternion<Scalar>& Quaternion<Scalar>::setFromTwoVectors(const MatrixBas { Vector3 v0 = a.normalized(); Vector3 v1 = b.normalized(); - Vector3 axis = v0.cross(v1); Scalar c = v0.dot(v1); // if dot == 1, vectors are the same @@ -353,7 +352,17 @@ inline Quaternion<Scalar>& Quaternion<Scalar>::setFromTwoVectors(const MatrixBas { // set to identity this->w() = 1; this->vec().setZero(); + return *this; } + // if dot == -1, vectors are opposites + if (ei_isApprox(c,Scalar(-1))) + { + this->vec() = v0.unitOrthogonal(); + this->w() = 0; + return *this; + } + + Vector3 axis = v0.cross(v1); Scalar s = ei_sqrt((Scalar(1)+c)*Scalar(2)); Scalar invs = Scalar(1)/s; this->vec() = axis * invs; |