Geometry_SSE.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009 Rohit Garg <rpg.314@gmail.com>
5 // Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
6 //
7 // Eigen is free software; you can redistribute it and/or
8 // modify it under the terms of the GNU Lesser General Public
9 // License as published by the Free Software Foundation; either
10 // version 3 of the License, or (at your option) any later version.
11 //
12 // Alternatively, you can redistribute it and/or
13 // modify it under the terms of the GNU General Public License as
14 // published by the Free Software Foundation; either version 2 of
15 // the License, or (at your option) any later version.
16 //
17 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
18 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
20 // GNU General Public License for more details.
21 //
22 // You should have received a copy of the GNU Lesser General Public
23 // License and a copy of the GNU General Public License along with
24 // Eigen. If not, see <http://www.gnu.org/licenses/>.
25 
26 #ifndef EIGEN_GEOMETRY_SSE_H
27 #define EIGEN_GEOMETRY_SSE_H
28 
29 namespace Eigen {
30 
31 namespace internal {
32 
33 template<class Derived, class OtherDerived>
34 struct quat_product<Architecture::SSE, Derived, OtherDerived, float, Aligned>
35 {
36  static inline Quaternion<float> run(const QuaternionBase<Derived>& _a, const QuaternionBase<OtherDerived>& _b)
37  {
38  const __m128 mask = _mm_castsi128_ps(_mm_setr_epi32(0,0,0,0x80000000));
39  Quaternion<float> res;
40  __m128 a = _a.coeffs().template packet<Aligned>(0);
41  __m128 b = _b.coeffs().template packet<Aligned>(0);
42  __m128 flip1 = _mm_xor_ps(_mm_mul_ps(vec4f_swizzle1(a,1,2,0,2),
43  vec4f_swizzle1(b,2,0,1,2)),mask);
44  __m128 flip2 = _mm_xor_ps(_mm_mul_ps(vec4f_swizzle1(a,3,3,3,1),
45  vec4f_swizzle1(b,0,1,2,1)),mask);
46  pstore(&res.x(),
47  _mm_add_ps(_mm_sub_ps(_mm_mul_ps(a,vec4f_swizzle1(b,3,3,3,3)),
48  _mm_mul_ps(vec4f_swizzle1(a,2,0,1,0),
49  vec4f_swizzle1(b,1,2,0,0))),
50  _mm_add_ps(flip1,flip2)));
51  return res;
52  }
53 };
54 
55 template<typename VectorLhs,typename VectorRhs>
56 struct cross3_impl<Architecture::SSE,VectorLhs,VectorRhs,float,true>
57 {
58  static inline typename plain_matrix_type<VectorLhs>::type
59  run(const VectorLhs& lhs, const VectorRhs& rhs)
60  {
61  __m128 a = lhs.template packet<VectorLhs::Flags&AlignedBit ? Aligned : Unaligned>(0);
62  __m128 b = rhs.template packet<VectorRhs::Flags&AlignedBit ? Aligned : Unaligned>(0);
63  __m128 mul1=_mm_mul_ps(vec4f_swizzle1(a,1,2,0,3),vec4f_swizzle1(b,2,0,1,3));
64  __m128 mul2=_mm_mul_ps(vec4f_swizzle1(a,2,0,1,3),vec4f_swizzle1(b,1,2,0,3));
65  typename plain_matrix_type<VectorLhs>::type res;
66  pstore(&res.x(),_mm_sub_ps(mul1,mul2));
67  return res;
68  }
69 };
70 
71 
72 
73 
74 template<class Derived, class OtherDerived>
75 struct quat_product<Architecture::SSE, Derived, OtherDerived, double, Aligned>
76 {
77  static inline Quaternion<double> run(const QuaternionBase<Derived>& _a, const QuaternionBase<OtherDerived>& _b)
78  {
79  const Packet2d mask = _mm_castsi128_pd(_mm_set_epi32(0x0,0x0,0x80000000,0x0));
80 
81  Quaternion<double> res;
82 
83  const double* a = _a.coeffs().data();
84  Packet2d b_xy = _b.coeffs().template packet<Aligned>(0);
85  Packet2d b_zw = _b.coeffs().template packet<Aligned>(2);
86  Packet2d a_xx = pset1<Packet2d>(a[0]);
87  Packet2d a_yy = pset1<Packet2d>(a[1]);
88  Packet2d a_zz = pset1<Packet2d>(a[2]);
89  Packet2d a_ww = pset1<Packet2d>(a[3]);
90 
91  // two temporaries:
92  Packet2d t1, t2;
93 
94  /*
95  * t1 = ww*xy + yy*zw
96  * t2 = zz*xy - xx*zw
97  * res.xy = t1 +/- swap(t2)
98  */
99  t1 = padd(pmul(a_ww, b_xy), pmul(a_yy, b_zw));
100  t2 = psub(pmul(a_zz, b_xy), pmul(a_xx, b_zw));
101 #ifdef EIGEN_VECTORIZE_SSE3
103  pstore(&res.x(), _mm_addsub_pd(t1, preverse(t2)));
104 #else
105  pstore(&res.x(), padd(t1, pxor(mask,preverse(t2))));
106 #endif
107 
108  /*
109  * t1 = ww*zw - yy*xy
110  * t2 = zz*zw + xx*xy
111  * res.zw = t1 -/+ swap(t2) = swap( swap(t1) +/- t2)
112  */
113  t1 = psub(pmul(a_ww, b_zw), pmul(a_yy, b_xy));
114  t2 = padd(pmul(a_zz, b_zw), pmul(a_xx, b_xy));
115 #ifdef EIGEN_VECTORIZE_SSE3
117  pstore(&res.z(), preverse(_mm_addsub_pd(preverse(t1), t2)));
118 #else
119  pstore(&res.z(), psub(t1, pxor(mask,preverse(t2))));
120 #endif
121 
122  return res;
123 }
124 };
125 
126 } // end namespace internal
127 
128 } // end namespace Eigen
129 
130 #endif // EIGEN_GEOMETRY_SSE_H