Eigen  3.2.91
Geometry_SSE.h
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 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_GEOMETRY_SSE_H
12 #define EIGEN_GEOMETRY_SSE_H
13 
14 namespace Eigen {
15 
16 namespace internal {
17 
18 template<class Derived, class OtherDerived>
19 struct quat_product<Architecture::SSE, Derived, OtherDerived, float, Aligned16>
20 {
21  static inline Quaternion<float> run(const QuaternionBase<Derived>& _a, const QuaternionBase<OtherDerived>& _b)
22  {
23  Quaternion<float> res;
24  const __m128 mask = _mm_setr_ps(0.f,0.f,0.f,-0.f);
25  __m128 a = _a.coeffs().template packet<Aligned16>(0);
26  __m128 b = _b.coeffs().template packet<Aligned16>(0);
27  __m128 s1 = _mm_mul_ps(vec4f_swizzle1(a,1,2,0,2),vec4f_swizzle1(b,2,0,1,2));
28  __m128 s2 = _mm_mul_ps(vec4f_swizzle1(a,3,3,3,1),vec4f_swizzle1(b,0,1,2,1));
29  pstore(&res.x(),
30  _mm_add_ps(_mm_sub_ps(_mm_mul_ps(a,vec4f_swizzle1(b,3,3,3,3)),
31  _mm_mul_ps(vec4f_swizzle1(a,2,0,1,0),
32  vec4f_swizzle1(b,1,2,0,0))),
33  _mm_xor_ps(mask,_mm_add_ps(s1,s2))));
34 
35  return res;
36  }
37 };
38 
39 template<class Derived, int Alignment>
40 struct quat_conj<Architecture::SSE, Derived, float, Alignment>
41 {
42  static inline Quaternion<float> run(const QuaternionBase<Derived>& q)
43  {
44  Quaternion<float> res;
45  const __m128 mask = _mm_setr_ps(-0.f,-0.f,-0.f,0.f);
46  pstore(&res.x(), _mm_xor_ps(mask, q.coeffs().template packet<Alignment>(0)));
47  return res;
48  }
49 };
50 
51 
52 template<typename VectorLhs,typename VectorRhs>
53 struct cross3_impl<Architecture::SSE,VectorLhs,VectorRhs,float,true>
54 {
55  static inline typename plain_matrix_type<VectorLhs>::type
56  run(const VectorLhs& lhs, const VectorRhs& rhs)
57  {
58  __m128 a = lhs.template packet<traits<VectorLhs>::Alignment>(0);
59  __m128 b = rhs.template packet<traits<VectorRhs>::Alignment>(0);
60  __m128 mul1=_mm_mul_ps(vec4f_swizzle1(a,1,2,0,3),vec4f_swizzle1(b,2,0,1,3));
61  __m128 mul2=_mm_mul_ps(vec4f_swizzle1(a,2,0,1,3),vec4f_swizzle1(b,1,2,0,3));
62  typename plain_matrix_type<VectorLhs>::type res;
63  pstore(&res.x(),_mm_sub_ps(mul1,mul2));
64  return res;
65  }
66 };
67 
68 
69 
70 
71 template<class Derived, class OtherDerived, int Alignment>
72 struct quat_product<Architecture::SSE, Derived, OtherDerived, double, Alignment>
73 {
74  static inline Quaternion<double> run(const QuaternionBase<Derived>& _a, const QuaternionBase<OtherDerived>& _b)
75  {
76  const Packet2d mask = _mm_castsi128_pd(_mm_set_epi32(0x0,0x0,0x80000000,0x0));
77 
78  Quaternion<double> res;
79 
80  const double* a = _a.coeffs().data();
81  Packet2d b_xy = _b.coeffs().template packet<Alignment>(0);
82  Packet2d b_zw = _b.coeffs().template packet<Alignment>(2);
83  Packet2d a_xx = pset1<Packet2d>(a[0]);
84  Packet2d a_yy = pset1<Packet2d>(a[1]);
85  Packet2d a_zz = pset1<Packet2d>(a[2]);
86  Packet2d a_ww = pset1<Packet2d>(a[3]);
87 
88  // two temporaries:
89  Packet2d t1, t2;
90 
91  /*
92  * t1 = ww*xy + yy*zw
93  * t2 = zz*xy - xx*zw
94  * res.xy = t1 +/- swap(t2)
95  */
96  t1 = padd(pmul(a_ww, b_xy), pmul(a_yy, b_zw));
97  t2 = psub(pmul(a_zz, b_xy), pmul(a_xx, b_zw));
98 #ifdef EIGEN_VECTORIZE_SSE3
99  EIGEN_UNUSED_VARIABLE(mask)
100  pstore(&res.x(), _mm_addsub_pd(t1, preverse(t2)));
101 #else
102  pstore(&res.x(), padd(t1, pxor(mask,preverse(t2))));
103 #endif
104 
105  /*
106  * t1 = ww*zw - yy*xy
107  * t2 = zz*zw + xx*xy
108  * res.zw = t1 -/+ swap(t2) = swap( swap(t1) +/- t2)
109  */
110  t1 = psub(pmul(a_ww, b_zw), pmul(a_yy, b_xy));
111  t2 = padd(pmul(a_zz, b_zw), pmul(a_xx, b_xy));
112 #ifdef EIGEN_VECTORIZE_SSE3
113  EIGEN_UNUSED_VARIABLE(mask)
114  pstore(&res.z(), preverse(_mm_addsub_pd(preverse(t1), t2)));
115 #else
116  pstore(&res.z(), psub(t1, pxor(mask,preverse(t2))));
117 #endif
118 
119  return res;
120 }
121 };
122 
123 template<class Derived, int Alignment>
124 struct quat_conj<Architecture::SSE, Derived, double, Alignment>
125 {
126  static inline Quaternion<double> run(const QuaternionBase<Derived>& q)
127  {
128  Quaternion<double> res;
129  const __m128d mask0 = _mm_setr_pd(-0.,-0.);
130  const __m128d mask2 = _mm_setr_pd(-0.,0.);
131  pstore(&res.x(), _mm_xor_pd(mask0, q.coeffs().template packet<Alignment>(0)));
132  pstore(&res.z(), _mm_xor_pd(mask2, q.coeffs().template packet<Alignment>(2)));
133  return res;
134  }
135 };
136 
137 } // end namespace internal
138 
139 } // end namespace Eigen
140 
141 #endif // EIGEN_GEOMETRY_SSE_H
Definition: LDLT.h:16
Definition: Constants.h:222
Definition: Eigen_Colamd.h:54