1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
|
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2009 Gael Guennebaud <g.gael@free.fr>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 3 of the License, or (at your option) any later version.
//
// Alternatively, you can redistribute it and/or
// modify it under the terms of the GNU General Public License as
// published by the Free Software Foundation; either version 2 of
// the License, or (at your option) any later version.
//
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.
#ifndef EIGEN_ALIGNED_VECTOR3
#define EIGEN_ALIGNED_VECTOR3
#include <Eigen/Geometry>
namespace Eigen {
/** \ingroup Unsupported_modules
* \defgroup AlignedVector3_Module Aligned vector3 module
*
* \code
* #include <unsupported/Eigen/AlignedVector3>
* \endcode
*/
//@{
/** \class AlignedVector3
*
* \brief A vectorization friendly 3D vector
*
* This class represents a 3D vector internally using a 4D vector
* such that vectorization can be seamlessly enabled. Of course,
* the same result can be achieved by directly using a 4D vector.
* This class makes this process simpler.
*
*/
// TODO specialize Cwise
template<typename _Scalar> class AlignedVector3;
template<typename _Scalar> struct ei_traits<AlignedVector3<_Scalar> >
: ei_traits<Matrix<_Scalar,3,1,0,4,1> >
{
};
template<typename _Scalar> class AlignedVector3
: public MatrixBase<AlignedVector3<_Scalar> >
{
typedef Matrix<_Scalar,4,1> CoeffType;
CoeffType m_coeffs;
public:
EIGEN_GENERIC_PUBLIC_INTERFACE(AlignedVector3)
using Base::operator*;
inline int rows() const { return 3; }
inline int cols() const { return 1; }
inline const Scalar& coeff(int row, int col) const
{ return m_coeffs.coeff(row, col); }
inline Scalar& coeffRef(int row, int col)
{ return m_coeffs.coeffRef(row, col); }
inline const Scalar& coeff(int index) const
{ return m_coeffs.coeff(index); }
inline Scalar& coeffRef(int index)
{ return m_coeffs.coeffRef(index);}
inline AlignedVector3(const Scalar& x, const Scalar& y, const Scalar& z)
: m_coeffs(x, y, z, Scalar(0))
{}
inline AlignedVector3(const AlignedVector3& other)
: Base(), m_coeffs(other.m_coeffs)
{}
template<typename XprType, int Size=XprType::SizeAtCompileTime>
struct generic_assign_selector {};
template<typename XprType> struct generic_assign_selector<XprType,4>
{
inline static void run(AlignedVector3& dest, const XprType& src)
{
dest.m_coeffs = src;
}
};
template<typename XprType> struct generic_assign_selector<XprType,3>
{
inline static void run(AlignedVector3& dest, const XprType& src)
{
dest.m_coeffs.template start<3>() = src;
dest.m_coeffs.w() = Scalar(0);
}
};
template<typename Derived>
inline explicit AlignedVector3(const MatrixBase<Derived>& other)
{
generic_assign_selector<Derived>::run(*this,other.derived());
}
inline AlignedVector3& operator=(const AlignedVector3& other)
{ m_coeffs = other.m_coeffs; return *this; }
inline AlignedVector3 operator+(const AlignedVector3& other) const
{ return AlignedVector3(m_coeffs + other.m_coeffs); }
inline AlignedVector3& operator+=(const AlignedVector3& other)
{ m_coeffs += other.m_coeffs; return *this; }
inline AlignedVector3 operator-(const AlignedVector3& other) const
{ return AlignedVector3(m_coeffs - other.m_coeffs); }
inline AlignedVector3 operator-=(const AlignedVector3& other)
{ m_coeffs -= other.m_coeffs; return *this; }
inline AlignedVector3 operator*(const Scalar& s) const
{ return AlignedVector3(m_coeffs * s); }
inline friend AlignedVector3 operator*(const Scalar& s,const AlignedVector3& vec)
{ return AlignedVector3(s * vec.m_coeffs); }
inline AlignedVector3& operator*=(const Scalar& s)
{ m_coeffs *= s; return *this; }
inline AlignedVector3 operator/(const Scalar& s) const
{ return AlignedVector3(m_coeffs / s); }
inline AlignedVector3& operator/=(const Scalar& s)
{ m_coeffs /= s; return *this; }
inline Scalar dot(const AlignedVector3& other) const
{
ei_assert(m_coeffs.w()==Scalar(0));
ei_assert(other.m_coeffs.w()==Scalar(0));
return m_coeffs.dot(other.m_coeffs);
}
inline void normalize()
{
m_coeffs /= norm();
}
inline AlignedVector3 normalized()
{
return AlignedVector3(m_coeffs / norm());
}
inline Scalar sum() const
{
ei_assert(m_coeffs.w()==Scalar(0));
return m_coeffs.sum();
}
inline Scalar squaredNorm() const
{
ei_assert(m_coeffs.w()==Scalar(0));
return m_coeffs.squaredNorm();
}
inline Scalar norm() const
{
return ei_sqrt(squaredNorm());
}
inline AlignedVector3 cross(const AlignedVector3& other) const
{
return AlignedVector3(m_coeffs.cross3(other.m_coeffs));
}
template<typename Derived>
inline bool isApprox(const MatrixBase<Derived>& other, RealScalar eps=precision<Scalar>()) const
{
return m_coeffs.template start<3>().isApprox(other,eps);
}
};
//@}
}
#endif // EIGEN_ALIGNED_VECTOR3
|