![]() |
Eigen
3.3.3
|
00001 // This file is part of Eigen, a lightweight C++ template library 00002 // for linear algebra. 00003 // 00004 // Copyright (C) 2009 Rohit Garg <rpg.314@gmail.com> 00005 // Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr> 00006 // 00007 // This Source Code Form is subject to the terms of the Mozilla 00008 // Public License v. 2.0. If a copy of the MPL was not distributed 00009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 00010 00011 #ifndef EIGEN_GEOMETRY_SSE_H 00012 #define EIGEN_GEOMETRY_SSE_H 00013 00014 namespace Eigen { 00015 00016 namespace internal { 00017 00018 template<class Derived, class OtherDerived> 00019 struct quat_product<Architecture::SSE, Derived, OtherDerived, float, Aligned16> 00020 { 00021 static inline Quaternion<float> run(const QuaternionBase<Derived>& _a, const QuaternionBase<OtherDerived>& _b) 00022 { 00023 Quaternion<float> res; 00024 const __m128 mask = _mm_setr_ps(0.f,0.f,0.f,-0.f); 00025 __m128 a = _a.coeffs().template packet<Aligned16>(0); 00026 __m128 b = _b.coeffs().template packet<Aligned16>(0); 00027 __m128 s1 = _mm_mul_ps(vec4f_swizzle1(a,1,2,0,2),vec4f_swizzle1(b,2,0,1,2)); 00028 __m128 s2 = _mm_mul_ps(vec4f_swizzle1(a,3,3,3,1),vec4f_swizzle1(b,0,1,2,1)); 00029 pstore(&res.x(), 00030 _mm_add_ps(_mm_sub_ps(_mm_mul_ps(a,vec4f_swizzle1(b,3,3,3,3)), 00031 _mm_mul_ps(vec4f_swizzle1(a,2,0,1,0), 00032 vec4f_swizzle1(b,1,2,0,0))), 00033 _mm_xor_ps(mask,_mm_add_ps(s1,s2)))); 00034 00035 return res; 00036 } 00037 }; 00038 00039 template<class Derived, int Alignment> 00040 struct quat_conj<Architecture::SSE, Derived, float, Alignment> 00041 { 00042 static inline Quaternion<float> run(const QuaternionBase<Derived>& q) 00043 { 00044 Quaternion<float> res; 00045 const __m128 mask = _mm_setr_ps(-0.f,-0.f,-0.f,0.f); 00046 pstore(&res.x(), _mm_xor_ps(mask, q.coeffs().template packet<Alignment>(0))); 00047 return res; 00048 } 00049 }; 00050 00051 00052 template<typename VectorLhs,typename VectorRhs> 00053 struct cross3_impl<Architecture::SSE,VectorLhs,VectorRhs,float,true> 00054 { 00055 static inline typename plain_matrix_type<VectorLhs>::type 00056 run(const VectorLhs& lhs, const VectorRhs& rhs) 00057 { 00058 __m128 a = lhs.template packet<traits<VectorLhs>::Alignment>(0); 00059 __m128 b = rhs.template packet<traits<VectorRhs>::Alignment>(0); 00060 __m128 mul1=_mm_mul_ps(vec4f_swizzle1(a,1,2,0,3),vec4f_swizzle1(b,2,0,1,3)); 00061 __m128 mul2=_mm_mul_ps(vec4f_swizzle1(a,2,0,1,3),vec4f_swizzle1(b,1,2,0,3)); 00062 typename plain_matrix_type<VectorLhs>::type res; 00063 pstore(&res.x(),_mm_sub_ps(mul1,mul2)); 00064 return res; 00065 } 00066 }; 00067 00068 00069 00070 00071 template<class Derived, class OtherDerived, int Alignment> 00072 struct quat_product<Architecture::SSE, Derived, OtherDerived, double, Alignment> 00073 { 00074 static inline Quaternion<double> run(const QuaternionBase<Derived>& _a, const QuaternionBase<OtherDerived>& _b) 00075 { 00076 const Packet2d mask = _mm_castsi128_pd(_mm_set_epi32(0x0,0x0,0x80000000,0x0)); 00077 00078 Quaternion<double> res; 00079 00080 const double* a = _a.coeffs().data(); 00081 Packet2d b_xy = _b.coeffs().template packet<Alignment>(0); 00082 Packet2d b_zw = _b.coeffs().template packet<Alignment>(2); 00083 Packet2d a_xx = pset1<Packet2d>(a[0]); 00084 Packet2d a_yy = pset1<Packet2d>(a[1]); 00085 Packet2d a_zz = pset1<Packet2d>(a[2]); 00086 Packet2d a_ww = pset1<Packet2d>(a[3]); 00087 00088 // two temporaries: 00089 Packet2d t1, t2; 00090 00091 /* 00092 * t1 = ww*xy + yy*zw 00093 * t2 = zz*xy - xx*zw 00094 * res.xy = t1 +/- swap(t2) 00095 */ 00096 t1 = padd(pmul(a_ww, b_xy), pmul(a_yy, b_zw)); 00097 t2 = psub(pmul(a_zz, b_xy), pmul(a_xx, b_zw)); 00098 #ifdef EIGEN_VECTORIZE_SSE3 00099 EIGEN_UNUSED_VARIABLE(mask) 00100 pstore(&res.x(), _mm_addsub_pd(t1, preverse(t2))); 00101 #else 00102 pstore(&res.x(), padd(t1, pxor(mask,preverse(t2)))); 00103 #endif 00104 00105 /* 00106 * t1 = ww*zw - yy*xy 00107 * t2 = zz*zw + xx*xy 00108 * res.zw = t1 -/+ swap(t2) = swap( swap(t1) +/- t2) 00109 */ 00110 t1 = psub(pmul(a_ww, b_zw), pmul(a_yy, b_xy)); 00111 t2 = padd(pmul(a_zz, b_zw), pmul(a_xx, b_xy)); 00112 #ifdef EIGEN_VECTORIZE_SSE3 00113 EIGEN_UNUSED_VARIABLE(mask) 00114 pstore(&res.z(), preverse(_mm_addsub_pd(preverse(t1), t2))); 00115 #else 00116 pstore(&res.z(), psub(t1, pxor(mask,preverse(t2)))); 00117 #endif 00118 00119 return res; 00120 } 00121 }; 00122 00123 template<class Derived, int Alignment> 00124 struct quat_conj<Architecture::SSE, Derived, double, Alignment> 00125 { 00126 static inline Quaternion<double> run(const QuaternionBase<Derived>& q) 00127 { 00128 Quaternion<double> res; 00129 const __m128d mask0 = _mm_setr_pd(-0.,-0.); 00130 const __m128d mask2 = _mm_setr_pd(-0.,0.); 00131 pstore(&res.x(), _mm_xor_pd(mask0, q.coeffs().template packet<Alignment>(0))); 00132 pstore(&res.z(), _mm_xor_pd(mask2, q.coeffs().template packet<Alignment>(2))); 00133 return res; 00134 } 00135 }; 00136 00137 } // end namespace internal 00138 00139 } // end namespace Eigen 00140 00141 #endif // EIGEN_GEOMETRY_SSE_H