11 #ifndef EIGEN_GEOMETRY_SIMD_H 12 #define EIGEN_GEOMETRY_SIMD_H 18 template<
class Derived,
class OtherDerived>
19 struct quat_product<Architecture::Target, Derived, OtherDerived, float>
22 AAlignment = traits<Derived>::Alignment,
23 BAlignment = traits<OtherDerived>::Alignment,
24 ResAlignment = traits<Quaternion<float> >::Alignment
26 static inline Quaternion<float> run(
const QuaternionBase<Derived>& _a,
const QuaternionBase<OtherDerived>& _b)
28 evaluator<typename Derived::Coefficients> ae(_a.coeffs());
29 evaluator<typename OtherDerived::Coefficients> be(_b.coeffs());
30 Quaternion<float> res;
31 const float neg_zero = numext::bit_cast<
float>(0x80000000u);
32 const float arr[4] = {0.f, 0.f, 0.f, neg_zero};
33 const Packet4f mask = ploadu<Packet4f>(arr);
34 Packet4f a = ae.template packet<AAlignment,Packet4f>(0);
35 Packet4f b = be.template packet<BAlignment,Packet4f>(0);
36 Packet4f s1 = pmul(vec4f_swizzle1(a,1,2,0,2),vec4f_swizzle1(b,2,0,1,2));
37 Packet4f s2 = pmul(vec4f_swizzle1(a,3,3,3,1),vec4f_swizzle1(b,0,1,2,1));
38 pstoret<float,Packet4f,ResAlignment>(
40 padd(psub(pmul(a,vec4f_swizzle1(b,3,3,3,3)),
41 pmul(vec4f_swizzle1(a,2,0,1,0),
42 vec4f_swizzle1(b,1,2,0,0))),
43 pxor(mask,padd(s1,s2))));
49 template<
class Derived>
50 struct quat_conj<Architecture::Target, Derived, float>
53 ResAlignment = traits<Quaternion<float> >::Alignment
55 static inline Quaternion<float> run(
const QuaternionBase<Derived>& q)
57 evaluator<typename Derived::Coefficients> qe(q.coeffs());
58 Quaternion<float> res;
59 const float neg_zero = numext::bit_cast<
float>(0x80000000u);
60 const float arr[4] = {neg_zero, neg_zero, neg_zero,0.f};
61 const Packet4f mask = ploadu<Packet4f>(arr);
62 pstoret<float,Packet4f,ResAlignment>(&res.x(), pxor(mask, qe.template packet<traits<Derived>::Alignment,Packet4f>(0)));
68 template<
typename VectorLhs,
typename VectorRhs>
69 struct cross3_impl<Architecture::Target,VectorLhs,VectorRhs,float,true>
72 ResAlignment = traits<typename plain_matrix_type<VectorLhs>::type>::Alignment
74 static inline typename plain_matrix_type<VectorLhs>::type
75 run(
const VectorLhs& lhs,
const VectorRhs& rhs)
77 evaluator<VectorLhs> lhs_eval(lhs);
78 evaluator<VectorRhs> rhs_eval(rhs);
79 Packet4f a = lhs_eval.template packet<traits<VectorLhs>::Alignment,Packet4f>(0);
80 Packet4f b = rhs_eval.template packet<traits<VectorRhs>::Alignment,Packet4f>(0);
81 Packet4f mul1 = pmul(vec4f_swizzle1(a,1,2,0,3),vec4f_swizzle1(b,2,0,1,3));
82 Packet4f mul2 = pmul(vec4f_swizzle1(a,2,0,1,3),vec4f_swizzle1(b,1,2,0,3));
83 typename plain_matrix_type<VectorLhs>::type res;
84 pstoret<float,Packet4f,ResAlignment>(&res.x(),psub(mul1,mul2));
91 #if (defined EIGEN_VECTORIZE_SSE) || (EIGEN_ARCH_ARM64) 93 template<
class Derived,
class OtherDerived>
94 struct quat_product<Architecture::Target, Derived, OtherDerived, double>
97 BAlignment = traits<OtherDerived>::Alignment,
98 ResAlignment = traits<Quaternion<double> >::Alignment
101 static inline Quaternion<double> run(
const QuaternionBase<Derived>& _a,
const QuaternionBase<OtherDerived>& _b)
103 Quaternion<double> res;
105 evaluator<typename Derived::Coefficients> ae(_a.coeffs());
106 evaluator<typename OtherDerived::Coefficients> be(_b.coeffs());
108 const double* a = _a.coeffs().data();
109 Packet2d b_xy = be.template packet<BAlignment,Packet2d>(0);
110 Packet2d b_zw = be.template packet<BAlignment,Packet2d>(2);
111 Packet2d a_xx = pset1<Packet2d>(a[0]);
112 Packet2d a_yy = pset1<Packet2d>(a[1]);
113 Packet2d a_zz = pset1<Packet2d>(a[2]);
114 Packet2d a_ww = pset1<Packet2d>(a[3]);
124 t1 = padd(pmul(a_ww, b_xy), pmul(a_yy, b_zw));
125 t2 = psub(pmul(a_zz, b_xy), pmul(a_xx, b_zw));
126 pstoret<double,Packet2d,ResAlignment>(&res.x(), paddsub(t1, preverse(t2)));
133 t1 = psub(pmul(a_ww, b_zw), pmul(a_yy, b_xy));
134 t2 = padd(pmul(a_zz, b_zw), pmul(a_xx, b_xy));
135 pstoret<double,Packet2d,ResAlignment>(&res.z(), preverse(paddsub(preverse(t1), t2)));
141 template<
class Derived>
142 struct quat_conj<Architecture::Target, Derived, double>
145 ResAlignment = traits<Quaternion<double> >::Alignment
147 static inline Quaternion<double> run(
const QuaternionBase<Derived>& q)
149 evaluator<typename Derived::Coefficients> qe(q.coeffs());
150 Quaternion<double> res;
151 const double neg_zero = numext::bit_cast<
double>(0x8000000000000000ull);
152 const double arr1[2] = {neg_zero, neg_zero};
153 const double arr2[2] = {neg_zero, 0.0};
154 const Packet2d mask0 = ploadu<Packet2d>(arr1);
155 const Packet2d mask2 = ploadu<Packet2d>(arr2);
156 pstoret<double,Packet2d,ResAlignment>(&res.x(), pxor(mask0, qe.template packet<traits<Derived>::Alignment,Packet2d>(0)));
157 pstoret<double,Packet2d,ResAlignment>(&res.z(), pxor(mask2, qe.template packet<traits<Derived>::Alignment,Packet2d>(2)));
162 #endif // end EIGEN_VECTORIZE_SSE_OR_EIGEN_ARCH_ARM64 168 #endif // EIGEN_GEOMETRY_SIMD_H Namespace containing all symbols from the Eigen library.
Definition: Core:141
Definition: Eigen_Colamd.h:50