1// This file is part of Eigen, a lightweight C++ template library
4// Copyright (C) 2009 Gael Guennebaud <g.gael@free.fr>
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/.
10#ifndef EIGEN_ALIGNED_VECTOR3
11#define EIGEN_ALIGNED_VECTOR3
13#include "../../Eigen/Geometry"
15#include "../../Eigen/src/Core/util/DisableStupidWarnings.h"
20 * \defgroup AlignedVector3_Module Aligned vector3 module
23 * #include <unsupported/Eigen/AlignedVector3>
29/** \class AlignedVector3
31 * \brief A vectorization friendly 3D vector
33 * This class represents a 3D vector internally using a 4D vector
34 * such that vectorization can be seamlessly enabled. Of course,
35 * the same result can be achieved by directly using a 4D vector.
36 * This class makes this process simpler.
39// TODO specialize Cwise
40template<typename _Scalar> class AlignedVector3;
43template<typename _Scalar> struct traits<AlignedVector3<_Scalar> >
44 : traits<Matrix<_Scalar,3,1,0,4,1> >
49template<typename _Scalar> class AlignedVector3
50 : public MatrixBase<AlignedVector3<_Scalar> >
52 typedef Matrix<_Scalar,4,1> CoeffType;
56 typedef MatrixBase<AlignedVector3<_Scalar> > Base;
57 EIGEN_DENSE_PUBLIC_INTERFACE(AlignedVector3)
58 using Base::operator*;
60 inline Index rows() const { return 3; }
61 inline Index cols() const { return 1; }
63 Scalar* data() { return m_coeffs.data(); }
64 const Scalar* data() const { return m_coeffs.data(); }
65 Index innerStride() const { return 1; }
66 Index outerStride() const { return 3; }
68 inline const Scalar& coeff(Index row, Index col) const
69 { return m_coeffs.coeff(row, col); }
71 inline Scalar& coeffRef(Index row, Index col)
72 { return m_coeffs.coeffRef(row, col); }
74 inline const Scalar& coeff(Index index) const
75 { return m_coeffs.coeff(index); }
77 inline Scalar& coeffRef(Index index)
78 { return m_coeffs.coeffRef(index);}
81 inline AlignedVector3()
84 inline AlignedVector3(const Scalar& x, const Scalar& y, const Scalar& z)
85 : m_coeffs(x, y, z, Scalar(0))
88 inline AlignedVector3(const AlignedVector3& other)
89 : Base(), m_coeffs(other.m_coeffs)
92 template<typename XprType, int Size=XprType::SizeAtCompileTime>
93 struct generic_assign_selector {};
95 template<typename XprType> struct generic_assign_selector<XprType,4>
97 inline static void run(AlignedVector3& dest, const XprType& src)
103 template<typename XprType> struct generic_assign_selector<XprType,3>
105 inline static void run(AlignedVector3& dest, const XprType& src)
107 dest.m_coeffs.template head<3>() = src;
108 dest.m_coeffs.w() = Scalar(0);
112 template<typename Derived>
113 inline AlignedVector3(const MatrixBase<Derived>& other)
115 generic_assign_selector<Derived>::run(*this,other.derived());
118 inline AlignedVector3& operator=(const AlignedVector3& other)
119 { m_coeffs = other.m_coeffs; return *this; }
121 template <typename Derived>
122 inline AlignedVector3& operator=(const MatrixBase<Derived>& other)
124 generic_assign_selector<Derived>::run(*this,other.derived());
128 inline AlignedVector3 operator+(const AlignedVector3& other) const
129 { return AlignedVector3(m_coeffs + other.m_coeffs); }
131 inline AlignedVector3& operator+=(const AlignedVector3& other)
132 { m_coeffs += other.m_coeffs; return *this; }
134 inline AlignedVector3 operator-(const AlignedVector3& other) const
135 { return AlignedVector3(m_coeffs - other.m_coeffs); }
137 inline AlignedVector3 operator-() const
138 { return AlignedVector3(-m_coeffs); }
140 inline AlignedVector3 operator-=(const AlignedVector3& other)
141 { m_coeffs -= other.m_coeffs; return *this; }
143 inline AlignedVector3 operator*(const Scalar& s) const
144 { return AlignedVector3(m_coeffs * s); }
146 inline friend AlignedVector3 operator*(const Scalar& s,const AlignedVector3& vec)
147 { return AlignedVector3(s * vec.m_coeffs); }
149 inline AlignedVector3& operator*=(const Scalar& s)
150 { m_coeffs *= s; return *this; }
152 inline AlignedVector3 operator/(const Scalar& s) const
153 { return AlignedVector3(m_coeffs / s); }
155 inline AlignedVector3& operator/=(const Scalar& s)
156 { m_coeffs /= s; return *this; }
158 inline Scalar dot(const AlignedVector3& other) const
160 eigen_assert(m_coeffs.w()==Scalar(0));
161 eigen_assert(other.m_coeffs.w()==Scalar(0));
162 return m_coeffs.dot(other.m_coeffs);
165 inline void normalize()
170 inline AlignedVector3 normalized() const
172 return AlignedVector3(m_coeffs / norm());
175 inline Scalar sum() const
177 eigen_assert(m_coeffs.w()==Scalar(0));
178 return m_coeffs.sum();
181 inline Scalar squaredNorm() const
183 eigen_assert(m_coeffs.w()==Scalar(0));
184 return m_coeffs.squaredNorm();
187 inline Scalar norm() const
190 return sqrt(squaredNorm());
193 inline AlignedVector3 cross(const AlignedVector3& other) const
195 return AlignedVector3(m_coeffs.cross3(other.m_coeffs));
198 template<typename Derived>
199 inline bool isApprox(const MatrixBase<Derived>& other, const RealScalar& eps=NumTraits<Scalar>::dummy_precision()) const
201 return m_coeffs.template head<3>().isApprox(other,eps);
204 CoeffType& coeffs() { return m_coeffs; }
205 const CoeffType& coeffs() const { return m_coeffs; }
210template<typename _Scalar>
211struct eval<AlignedVector3<_Scalar>, Dense>
213 typedef const AlignedVector3<_Scalar>& type;
216template<typename Scalar>
217struct evaluator<AlignedVector3<Scalar> >
218 : evaluator<Matrix<Scalar,4,1> >
220 typedef AlignedVector3<Scalar> XprType;
221 typedef evaluator<Matrix<Scalar,4,1> > Base;
223 evaluator(const XprType &m) : Base(m.coeffs()) {}
232#include "../../Eigen/src/Core/util/ReenableStupidWarnings.h"
234#endif // EIGEN_ALIGNED_VECTOR3