aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Geometry/OrthoMethods.h
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2009-02-17 09:53:05 +0000
committerGravatar Gael Guennebaud <g.gael@free.fr>2009-02-17 09:53:05 +0000
commite6f1104b57f19dff773b4f22d26d6aacabd1bdb2 (patch)
tree5de1c7a6d6eb47d894f119b5bda0107bd531000d /Eigen/src/Geometry/OrthoMethods.h
parent67b4fab4e30a59d9a7e001ef25938d1767371569 (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/OrthoMethods.h')
-rw-r--r--Eigen/src/Geometry/OrthoMethods.h38
1 files changed, 37 insertions, 1 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
{