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
10#ifndef EIGEN_ALIGNED_VECTOR3
11#define EIGEN_ALIGNED_VECTOR3
12
13#include <Eigen/Geometry>
14
15namespace Eigen {
16
17/**
18  * \defgroup AlignedVector3_Module Aligned vector3 module
19  *
20  * \code
21  * #include <unsupported/Eigen/AlignedVector3>
22  * \endcode
23  */
24  //@{
25
26
27/** \class AlignedVector3
28  *
29  * \brief A vectorization friendly 3D vector
30  *
31  * This class represents a 3D vector internally using a 4D vector
32  * such that vectorization can be seamlessly enabled. Of course,
33  * the same result can be achieved by directly using a 4D vector.
34  * This class makes this process simpler.
35  *
36  */
37// TODO specialize Cwise
38template<typename _Scalar> class AlignedVector3;
39
40namespace internal {
41template<typename _Scalar> struct traits<AlignedVector3<_Scalar> >
42  : traits<Matrix<_Scalar,3,1,0,4,1> >
43{
44};
45}
46
47template<typename _Scalar> class AlignedVector3
48  : public MatrixBase<AlignedVector3<_Scalar> >
49{
50    typedef Matrix<_Scalar,4,1> CoeffType;
51    CoeffType m_coeffs;
52  public:
53
54    typedef MatrixBase<AlignedVector3<_Scalar> > Base;
55    EIGEN_DENSE_PUBLIC_INTERFACE(AlignedVector3)
56    using Base::operator*;
57
58    inline Index rows() const { return 3; }
59    inline Index cols() const { return 1; }
60
61    Scalar* data() { return m_coeffs.data(); }
62    const Scalar* data() const { return m_coeffs.data(); }
63    Index innerStride() const { return 1; }
64    Index outerStride() const { return 3; }
65
66    inline const Scalar& coeff(Index row, Index col) const
67    { return m_coeffs.coeff(row, col); }
68
69    inline Scalar& coeffRef(Index row, Index col)
70    { return m_coeffs.coeffRef(row, col); }
71
72    inline const Scalar& coeff(Index index) const
73    { return m_coeffs.coeff(index); }
74
75    inline Scalar& coeffRef(Index index)
76    { return m_coeffs.coeffRef(index);}
77
78
79    inline AlignedVector3(const Scalar& x, const Scalar& y, const Scalar& z)
80      : m_coeffs(x, y, z, Scalar(0))
81    {}
82
83    inline AlignedVector3(const AlignedVector3& other)
84      : Base(), m_coeffs(other.m_coeffs)
85    {}
86
87    template<typename XprType, int Size=XprType::SizeAtCompileTime>
88    struct generic_assign_selector {};
89
90    template<typename XprType> struct generic_assign_selector<XprType,4>
91    {
92      inline static void run(AlignedVector3& dest, const XprType& src)
93      {
94        dest.m_coeffs = src;
95      }
96    };
97
98    template<typename XprType> struct generic_assign_selector<XprType,3>
99    {
100      inline static void run(AlignedVector3& dest, const XprType& src)
101      {
102        dest.m_coeffs.template head<3>() = src;
103        dest.m_coeffs.w() = Scalar(0);
104      }
105    };
106
107    template<typename Derived>
108    inline AlignedVector3(const MatrixBase<Derived>& other)
109    {
110      generic_assign_selector<Derived>::run(*this,other.derived());
111    }
112
113    inline AlignedVector3& operator=(const AlignedVector3& other)
114    { m_coeffs = other.m_coeffs; return *this; }
115
116    template <typename Derived>
117    inline AlignedVector3& operator=(const MatrixBase<Derived>& other)
118    {
119      generic_assign_selector<Derived>::run(*this,other.derived());
120      return *this;
121    }
122
123    inline AlignedVector3 operator+(const AlignedVector3& other) const
124    { return AlignedVector3(m_coeffs + other.m_coeffs); }
125
126    inline AlignedVector3& operator+=(const AlignedVector3& other)
127    { m_coeffs += other.m_coeffs; return *this; }
128
129    inline AlignedVector3 operator-(const AlignedVector3& other) const
130    { return AlignedVector3(m_coeffs - other.m_coeffs); }
131
132    inline AlignedVector3 operator-=(const AlignedVector3& other)
133    { m_coeffs -= other.m_coeffs; return *this; }
134
135    inline AlignedVector3 operator*(const Scalar& s) const
136    { return AlignedVector3(m_coeffs * s); }
137
138    inline friend AlignedVector3 operator*(const Scalar& s,const AlignedVector3& vec)
139    { return AlignedVector3(s * vec.m_coeffs); }
140
141    inline AlignedVector3& operator*=(const Scalar& s)
142    { m_coeffs *= s; return *this; }
143
144    inline AlignedVector3 operator/(const Scalar& s) const
145    { return AlignedVector3(m_coeffs / s); }
146
147    inline AlignedVector3& operator/=(const Scalar& s)
148    { m_coeffs /= s; return *this; }
149
150    inline Scalar dot(const AlignedVector3& other) const
151    {
152      eigen_assert(m_coeffs.w()==Scalar(0));
153      eigen_assert(other.m_coeffs.w()==Scalar(0));
154      return m_coeffs.dot(other.m_coeffs);
155    }
156
157    inline void normalize()
158    {
159      m_coeffs /= norm();
160    }
161
162    inline AlignedVector3 normalized() const
163    {
164      return AlignedVector3(m_coeffs / norm());
165    }
166
167    inline Scalar sum() const
168    {
169      eigen_assert(m_coeffs.w()==Scalar(0));
170      return m_coeffs.sum();
171    }
172
173    inline Scalar squaredNorm() const
174    {
175      eigen_assert(m_coeffs.w()==Scalar(0));
176      return m_coeffs.squaredNorm();
177    }
178
179    inline Scalar norm() const
180    {
181      using std::sqrt;
182      return sqrt(squaredNorm());
183    }
184
185    inline AlignedVector3 cross(const AlignedVector3& other) const
186    {
187      return AlignedVector3(m_coeffs.cross3(other.m_coeffs));
188    }
189
190    template<typename Derived>
191    inline bool isApprox(const MatrixBase<Derived>& other, const RealScalar& eps=NumTraits<Scalar>::dummy_precision()) const
192    {
193      return m_coeffs.template head<3>().isApprox(other,eps);
194    }
195
196    CoeffType& coeffs() { return m_coeffs; }
197    const CoeffType& coeffs() const { return m_coeffs; }
198};
199
200namespace internal {
201
202template<typename _Scalar>
203struct eval<AlignedVector3<_Scalar>, Dense>
204{
205 typedef const AlignedVector3<_Scalar>& type;
206};
207
208template<typename Scalar>
209struct evaluator<AlignedVector3<Scalar> >
210  : evaluator<Matrix<Scalar,4,1> >
211{
212  typedef AlignedVector3<Scalar> XprType;
213  typedef evaluator<Matrix<Scalar,4,1> > Base;
214
215  evaluator(const XprType &m) : Base(m.coeffs()) {}
216};
217
218}
219
220//@}
221
222}
223
224#endif // EIGEN_ALIGNED_VECTOR3
225