AlignedVector3 6.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2009 Gael Guennebaud <g.gael@free.fr>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #ifndef EIGEN_ALIGNED_VECTOR3
  10. #define EIGEN_ALIGNED_VECTOR3
  11. #include <Eigen/Geometry>
  12. namespace Eigen {
  13. /**
  14. * \defgroup AlignedVector3_Module Aligned vector3 module
  15. *
  16. * \code
  17. * #include <unsupported/Eigen/AlignedVector3>
  18. * \endcode
  19. */
  20. //@{
  21. /** \class AlignedVector3
  22. *
  23. * \brief A vectorization friendly 3D vector
  24. *
  25. * This class represents a 3D vector internally using a 4D vector
  26. * such that vectorization can be seamlessly enabled. Of course,
  27. * the same result can be achieved by directly using a 4D vector.
  28. * This class makes this process simpler.
  29. *
  30. */
  31. // TODO specialize Cwise
  32. template<typename _Scalar> class AlignedVector3;
  33. namespace internal {
  34. template<typename _Scalar> struct traits<AlignedVector3<_Scalar> >
  35. : traits<Matrix<_Scalar,3,1,0,4,1> >
  36. {
  37. };
  38. }
  39. template<typename _Scalar> class AlignedVector3
  40. : public MatrixBase<AlignedVector3<_Scalar> >
  41. {
  42. typedef Matrix<_Scalar,4,1> CoeffType;
  43. CoeffType m_coeffs;
  44. public:
  45. typedef MatrixBase<AlignedVector3<_Scalar> > Base;
  46. EIGEN_DENSE_PUBLIC_INTERFACE(AlignedVector3)
  47. using Base::operator*;
  48. inline Index rows() const { return 3; }
  49. inline Index cols() const { return 1; }
  50. Scalar* data() { return m_coeffs.data(); }
  51. const Scalar* data() const { return m_coeffs.data(); }
  52. Index innerStride() const { return 1; }
  53. Index outerStride() const { return 3; }
  54. inline const Scalar& coeff(Index row, Index col) const
  55. { return m_coeffs.coeff(row, col); }
  56. inline Scalar& coeffRef(Index row, Index col)
  57. { return m_coeffs.coeffRef(row, col); }
  58. inline const Scalar& coeff(Index index) const
  59. { return m_coeffs.coeff(index); }
  60. inline Scalar& coeffRef(Index index)
  61. { return m_coeffs.coeffRef(index);}
  62. inline AlignedVector3(const Scalar& x, const Scalar& y, const Scalar& z)
  63. : m_coeffs(x, y, z, Scalar(0))
  64. {}
  65. inline AlignedVector3(const AlignedVector3& other)
  66. : Base(), m_coeffs(other.m_coeffs)
  67. {}
  68. template<typename XprType, int Size=XprType::SizeAtCompileTime>
  69. struct generic_assign_selector {};
  70. template<typename XprType> struct generic_assign_selector<XprType,4>
  71. {
  72. inline static void run(AlignedVector3& dest, const XprType& src)
  73. {
  74. dest.m_coeffs = src;
  75. }
  76. };
  77. template<typename XprType> struct generic_assign_selector<XprType,3>
  78. {
  79. inline static void run(AlignedVector3& dest, const XprType& src)
  80. {
  81. dest.m_coeffs.template head<3>() = src;
  82. dest.m_coeffs.w() = Scalar(0);
  83. }
  84. };
  85. template<typename Derived>
  86. inline AlignedVector3(const MatrixBase<Derived>& other)
  87. {
  88. generic_assign_selector<Derived>::run(*this,other.derived());
  89. }
  90. inline AlignedVector3& operator=(const AlignedVector3& other)
  91. { m_coeffs = other.m_coeffs; return *this; }
  92. template <typename Derived>
  93. inline AlignedVector3& operator=(const MatrixBase<Derived>& other)
  94. {
  95. generic_assign_selector<Derived>::run(*this,other.derived());
  96. return *this;
  97. }
  98. inline AlignedVector3 operator+(const AlignedVector3& other) const
  99. { return AlignedVector3(m_coeffs + other.m_coeffs); }
  100. inline AlignedVector3& operator+=(const AlignedVector3& other)
  101. { m_coeffs += other.m_coeffs; return *this; }
  102. inline AlignedVector3 operator-(const AlignedVector3& other) const
  103. { return AlignedVector3(m_coeffs - other.m_coeffs); }
  104. inline AlignedVector3 operator-=(const AlignedVector3& other)
  105. { m_coeffs -= other.m_coeffs; return *this; }
  106. inline AlignedVector3 operator*(const Scalar& s) const
  107. { return AlignedVector3(m_coeffs * s); }
  108. inline friend AlignedVector3 operator*(const Scalar& s,const AlignedVector3& vec)
  109. { return AlignedVector3(s * vec.m_coeffs); }
  110. inline AlignedVector3& operator*=(const Scalar& s)
  111. { m_coeffs *= s; return *this; }
  112. inline AlignedVector3 operator/(const Scalar& s) const
  113. { return AlignedVector3(m_coeffs / s); }
  114. inline AlignedVector3& operator/=(const Scalar& s)
  115. { m_coeffs /= s; return *this; }
  116. inline Scalar dot(const AlignedVector3& other) const
  117. {
  118. eigen_assert(m_coeffs.w()==Scalar(0));
  119. eigen_assert(other.m_coeffs.w()==Scalar(0));
  120. return m_coeffs.dot(other.m_coeffs);
  121. }
  122. inline void normalize()
  123. {
  124. m_coeffs /= norm();
  125. }
  126. inline AlignedVector3 normalized() const
  127. {
  128. return AlignedVector3(m_coeffs / norm());
  129. }
  130. inline Scalar sum() const
  131. {
  132. eigen_assert(m_coeffs.w()==Scalar(0));
  133. return m_coeffs.sum();
  134. }
  135. inline Scalar squaredNorm() const
  136. {
  137. eigen_assert(m_coeffs.w()==Scalar(0));
  138. return m_coeffs.squaredNorm();
  139. }
  140. inline Scalar norm() const
  141. {
  142. using std::sqrt;
  143. return sqrt(squaredNorm());
  144. }
  145. inline AlignedVector3 cross(const AlignedVector3& other) const
  146. {
  147. return AlignedVector3(m_coeffs.cross3(other.m_coeffs));
  148. }
  149. template<typename Derived>
  150. inline bool isApprox(const MatrixBase<Derived>& other, const RealScalar& eps=NumTraits<Scalar>::dummy_precision()) const
  151. {
  152. return m_coeffs.template head<3>().isApprox(other,eps);
  153. }
  154. CoeffType& coeffs() { return m_coeffs; }
  155. const CoeffType& coeffs() const { return m_coeffs; }
  156. };
  157. namespace internal {
  158. template<typename _Scalar>
  159. struct eval<AlignedVector3<_Scalar>, Dense>
  160. {
  161. typedef const AlignedVector3<_Scalar>& type;
  162. };
  163. template<typename Scalar>
  164. struct evaluator<AlignedVector3<Scalar> >
  165. : evaluator<Matrix<Scalar,4,1> >
  166. {
  167. typedef AlignedVector3<Scalar> XprType;
  168. typedef evaluator<Matrix<Scalar,4,1> > Base;
  169. evaluator(const XprType &m) : Base(m.coeffs()) {}
  170. };
  171. }
  172. //@}
  173. }
  174. #endif // EIGEN_ALIGNED_VECTOR3