OrthoMethods.h 8.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
  6. //
  7. // This Source Code Form is subject to the terms of the Mozilla
  8. // Public License v. 2.0. If a copy of the MPL was not distributed
  9. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  10. #ifndef EIGEN_ORTHOMETHODS_H
  11. #define EIGEN_ORTHOMETHODS_H
  12. namespace Eigen {
  13. /** \geometry_module \ingroup Geometry_Module
  14. *
  15. * \returns the cross product of \c *this and \a other
  16. *
  17. * Here is a very good explanation of cross-product: http://xkcd.com/199/
  18. *
  19. * With complex numbers, the cross product is implemented as
  20. * \f$ (\mathbf{a}+i\mathbf{b}) \times (\mathbf{c}+i\mathbf{d}) = (\mathbf{a} \times \mathbf{c} - \mathbf{b} \times \mathbf{d}) - i(\mathbf{a} \times \mathbf{d} - \mathbf{b} \times \mathbf{c})\f$
  21. *
  22. * \sa MatrixBase::cross3()
  23. */
  24. template<typename Derived>
  25. template<typename OtherDerived>
  26. #ifndef EIGEN_PARSED_BY_DOXYGEN
  27. EIGEN_DEVICE_FUNC inline typename MatrixBase<Derived>::template cross_product_return_type<OtherDerived>::type
  28. #else
  29. inline typename MatrixBase<Derived>::PlainObject
  30. #endif
  31. MatrixBase<Derived>::cross(const MatrixBase<OtherDerived>& other) const
  32. {
  33. EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,3)
  34. EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3)
  35. // Note that there is no need for an expression here since the compiler
  36. // optimize such a small temporary very well (even within a complex expression)
  37. typename internal::nested_eval<Derived,2>::type lhs(derived());
  38. typename internal::nested_eval<OtherDerived,2>::type rhs(other.derived());
  39. return typename cross_product_return_type<OtherDerived>::type(
  40. numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),
  41. numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),
  42. numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0))
  43. );
  44. }
  45. namespace internal {
  46. template< int Arch,typename VectorLhs,typename VectorRhs,
  47. typename Scalar = typename VectorLhs::Scalar,
  48. bool Vectorizable = bool((VectorLhs::Flags&VectorRhs::Flags)&PacketAccessBit)>
  49. struct cross3_impl {
  50. EIGEN_DEVICE_FUNC static inline typename internal::plain_matrix_type<VectorLhs>::type
  51. run(const VectorLhs& lhs, const VectorRhs& rhs)
  52. {
  53. return typename internal::plain_matrix_type<VectorLhs>::type(
  54. numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),
  55. numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),
  56. numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0)),
  57. 0
  58. );
  59. }
  60. };
  61. }
  62. /** \geometry_module \ingroup Geometry_Module
  63. *
  64. * \returns the cross product of \c *this and \a other using only the x, y, and z coefficients
  65. *
  66. * The size of \c *this and \a other must be four. This function is especially useful
  67. * when using 4D vectors instead of 3D ones to get advantage of SSE/AltiVec vectorization.
  68. *
  69. * \sa MatrixBase::cross()
  70. */
  71. template<typename Derived>
  72. template<typename OtherDerived>
  73. EIGEN_DEVICE_FUNC inline typename MatrixBase<Derived>::PlainObject
  74. MatrixBase<Derived>::cross3(const MatrixBase<OtherDerived>& other) const
  75. {
  76. EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,4)
  77. EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,4)
  78. typedef typename internal::nested_eval<Derived,2>::type DerivedNested;
  79. typedef typename internal::nested_eval<OtherDerived,2>::type OtherDerivedNested;
  80. DerivedNested lhs(derived());
  81. OtherDerivedNested rhs(other.derived());
  82. return internal::cross3_impl<Architecture::Target,
  83. typename internal::remove_all<DerivedNested>::type,
  84. typename internal::remove_all<OtherDerivedNested>::type>::run(lhs,rhs);
  85. }
  86. /** \geometry_module \ingroup Geometry_Module
  87. *
  88. * \returns a matrix expression of the cross product of each column or row
  89. * of the referenced expression with the \a other vector.
  90. *
  91. * The referenced matrix must have one dimension equal to 3.
  92. * The result matrix has the same dimensions than the referenced one.
  93. *
  94. * \sa MatrixBase::cross() */
  95. template<typename ExpressionType, int Direction>
  96. template<typename OtherDerived>
  97. EIGEN_DEVICE_FUNC
  98. const typename VectorwiseOp<ExpressionType,Direction>::CrossReturnType
  99. VectorwiseOp<ExpressionType,Direction>::cross(const MatrixBase<OtherDerived>& other) const
  100. {
  101. EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3)
  102. EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename OtherDerived::Scalar>::value),
  103. YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
  104. typename internal::nested_eval<ExpressionType,2>::type mat(_expression());
  105. typename internal::nested_eval<OtherDerived,2>::type vec(other.derived());
  106. CrossReturnType res(_expression().rows(),_expression().cols());
  107. if(Direction==Vertical)
  108. {
  109. eigen_assert(CrossReturnType::RowsAtCompileTime==3 && "the matrix must have exactly 3 rows");
  110. res.row(0) = (mat.row(1) * vec.coeff(2) - mat.row(2) * vec.coeff(1)).conjugate();
  111. res.row(1) = (mat.row(2) * vec.coeff(0) - mat.row(0) * vec.coeff(2)).conjugate();
  112. res.row(2) = (mat.row(0) * vec.coeff(1) - mat.row(1) * vec.coeff(0)).conjugate();
  113. }
  114. else
  115. {
  116. eigen_assert(CrossReturnType::ColsAtCompileTime==3 && "the matrix must have exactly 3 columns");
  117. res.col(0) = (mat.col(1) * vec.coeff(2) - mat.col(2) * vec.coeff(1)).conjugate();
  118. res.col(1) = (mat.col(2) * vec.coeff(0) - mat.col(0) * vec.coeff(2)).conjugate();
  119. res.col(2) = (mat.col(0) * vec.coeff(1) - mat.col(1) * vec.coeff(0)).conjugate();
  120. }
  121. return res;
  122. }
  123. namespace internal {
  124. template<typename Derived, int Size = Derived::SizeAtCompileTime>
  125. struct unitOrthogonal_selector
  126. {
  127. typedef typename plain_matrix_type<Derived>::type VectorType;
  128. typedef typename traits<Derived>::Scalar Scalar;
  129. typedef typename NumTraits<Scalar>::Real RealScalar;
  130. typedef Matrix<Scalar,2,1> Vector2;
  131. EIGEN_DEVICE_FUNC
  132. static inline VectorType run(const Derived& src)
  133. {
  134. VectorType perp = VectorType::Zero(src.size());
  135. Index maxi = 0;
  136. Index sndi = 0;
  137. src.cwiseAbs().maxCoeff(&maxi);
  138. if (maxi==0)
  139. sndi = 1;
  140. RealScalar invnm = RealScalar(1)/(Vector2() << src.coeff(sndi),src.coeff(maxi)).finished().norm();
  141. perp.coeffRef(maxi) = -numext::conj(src.coeff(sndi)) * invnm;
  142. perp.coeffRef(sndi) = numext::conj(src.coeff(maxi)) * invnm;
  143. return perp;
  144. }
  145. };
  146. template<typename Derived>
  147. struct unitOrthogonal_selector<Derived,3>
  148. {
  149. typedef typename plain_matrix_type<Derived>::type VectorType;
  150. typedef typename traits<Derived>::Scalar Scalar;
  151. typedef typename NumTraits<Scalar>::Real RealScalar;
  152. EIGEN_DEVICE_FUNC
  153. static inline VectorType run(const Derived& src)
  154. {
  155. VectorType perp;
  156. /* Let us compute the crossed product of *this with a vector
  157. * that is not too close to being colinear to *this.
  158. */
  159. /* unless the x and y coords are both close to zero, we can
  160. * simply take ( -y, x, 0 ) and normalize it.
  161. */
  162. if((!isMuchSmallerThan(src.x(), src.z()))
  163. || (!isMuchSmallerThan(src.y(), src.z())))
  164. {
  165. RealScalar invnm = RealScalar(1)/src.template head<2>().norm();
  166. perp.coeffRef(0) = -numext::conj(src.y())*invnm;
  167. perp.coeffRef(1) = numext::conj(src.x())*invnm;
  168. perp.coeffRef(2) = 0;
  169. }
  170. /* if both x and y are close to zero, then the vector is close
  171. * to the z-axis, so it's far from colinear to the x-axis for instance.
  172. * So we take the crossed product with (1,0,0) and normalize it.
  173. */
  174. else
  175. {
  176. RealScalar invnm = RealScalar(1)/src.template tail<2>().norm();
  177. perp.coeffRef(0) = 0;
  178. perp.coeffRef(1) = -numext::conj(src.z())*invnm;
  179. perp.coeffRef(2) = numext::conj(src.y())*invnm;
  180. }
  181. return perp;
  182. }
  183. };
  184. template<typename Derived>
  185. struct unitOrthogonal_selector<Derived,2>
  186. {
  187. typedef typename plain_matrix_type<Derived>::type VectorType;
  188. EIGEN_DEVICE_FUNC
  189. static inline VectorType run(const Derived& src)
  190. { return VectorType(-numext::conj(src.y()), numext::conj(src.x())).normalized(); }
  191. };
  192. } // end namespace internal
  193. /** \geometry_module \ingroup Geometry_Module
  194. *
  195. * \returns a unit vector which is orthogonal to \c *this
  196. *
  197. * The size of \c *this must be at least 2. If the size is exactly 2,
  198. * then the returned vector is a counter clock wise rotation of \c *this, i.e., (-y,x).normalized().
  199. *
  200. * \sa cross()
  201. */
  202. template<typename Derived>
  203. EIGEN_DEVICE_FUNC typename MatrixBase<Derived>::PlainObject
  204. MatrixBase<Derived>::unitOrthogonal() const
  205. {
  206. EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
  207. return internal::unitOrthogonal_selector<Derived>::run(derived());
  208. }
  209. } // end namespace Eigen
  210. #endif // EIGEN_ORTHOMETHODS_H