$darkmode
Eigen  5.0.1-dev
Geometry_SIMD.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_SIMD_H
12 #define EIGEN_GEOMETRY_SIMD_H
13 
14 // IWYU pragma: private
15 #include "../InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
19 namespace internal {
20 
21 template <class Derived, class OtherDerived>
22 struct quat_product<Architecture::Target, Derived, OtherDerived, float> {
23  enum {
24  AAlignment = traits<Derived>::Alignment,
25  BAlignment = traits<OtherDerived>::Alignment,
26  ResAlignment = traits<Quaternion<float> >::Alignment
27  };
28  static inline Quaternion<float> run(const QuaternionBase<Derived>& _a, const QuaternionBase<OtherDerived>& _b) {
29  evaluator<typename Derived::Coefficients> ae(_a.coeffs());
30  evaluator<typename OtherDerived::Coefficients> be(_b.coeffs());
31  Quaternion<float> res;
32  const float neg_zero = numext::bit_cast<float>(0x80000000u);
33  const float arr[4] = {0.f, 0.f, 0.f, neg_zero};
34  const Packet4f mask = ploadu<Packet4f>(arr);
35  Packet4f a = ae.template packet<AAlignment, Packet4f>(0);
36  Packet4f b = be.template packet<BAlignment, Packet4f>(0);
37  Packet4f s1 = pmul(vec4f_swizzle1(a, 1, 2, 0, 2), vec4f_swizzle1(b, 2, 0, 1, 2));
38  Packet4f s2 = pmul(vec4f_swizzle1(a, 3, 3, 3, 1), vec4f_swizzle1(b, 0, 1, 2, 1));
39  pstoret<float, Packet4f, ResAlignment>(
40  &res.x(), padd(psub(pmul(a, vec4f_swizzle1(b, 3, 3, 3, 3)),
41  pmul(vec4f_swizzle1(a, 2, 0, 1, 0), vec4f_swizzle1(b, 1, 2, 0, 0))),
42  pxor(mask, padd(s1, s2))));
43 
44  return res;
45  }
46 };
47 
48 template <class Derived>
49 struct quat_conj<Architecture::Target, Derived, float> {
50  enum { ResAlignment = traits<Quaternion<float> >::Alignment };
51  static inline Quaternion<float> run(const QuaternionBase<Derived>& q) {
52  evaluator<typename Derived::Coefficients> qe(q.coeffs());
53  Quaternion<float> res;
54  const float neg_zero = numext::bit_cast<float>(0x80000000u);
55  const float arr[4] = {neg_zero, neg_zero, neg_zero, 0.f};
56  const Packet4f mask = ploadu<Packet4f>(arr);
57  pstoret<float, Packet4f, ResAlignment>(&res.x(),
58  pxor(mask, qe.template packet<traits<Derived>::Alignment, Packet4f>(0)));
59  return res;
60  }
61 };
62 
63 template <typename VectorLhs, typename VectorRhs>
64 struct cross3_impl<Architecture::Target, VectorLhs, VectorRhs, float, true> {
65  using DstPlainType = typename plain_matrix_type<VectorLhs>::type;
66  static constexpr int DstAlignment = evaluator<DstPlainType>::Alignment;
67  static constexpr int LhsAlignment = evaluator<VectorLhs>::Alignment;
68  static constexpr int RhsAlignment = evaluator<VectorRhs>::Alignment;
69  static inline DstPlainType run(const VectorLhs& lhs, const VectorRhs& rhs) {
70  evaluator<VectorLhs> lhs_eval(lhs);
71  evaluator<VectorRhs> rhs_eval(rhs);
72  Packet4f a = lhs_eval.template packet<LhsAlignment, Packet4f>(0);
73  Packet4f b = rhs_eval.template packet<RhsAlignment, Packet4f>(0);
74  Packet4f mul1 = pmul(vec4f_swizzle1(a, 1, 2, 0, 3), vec4f_swizzle1(b, 2, 0, 1, 3));
75  Packet4f mul2 = pmul(vec4f_swizzle1(a, 2, 0, 1, 3), vec4f_swizzle1(b, 1, 2, 0, 3));
76  DstPlainType res;
77  pstoret<float, Packet4f, DstAlignment>(res.data(), psub(mul1, mul2));
78  // Ensure last component is 0 in case original a or b contain inf/nan.
79  res[3] = 0.0f;
80  return res;
81  }
82 };
83 
84 #if (defined EIGEN_VECTORIZE_SSE) || (EIGEN_ARCH_ARM64)
85 
86 template <class Derived, class OtherDerived>
87 struct quat_product<Architecture::Target, Derived, OtherDerived, double> {
88  enum { BAlignment = traits<OtherDerived>::Alignment, ResAlignment = traits<Quaternion<double> >::Alignment };
89 
90  static inline Quaternion<double> run(const QuaternionBase<Derived>& _a, const QuaternionBase<OtherDerived>& _b) {
91  Quaternion<double> res;
92 
93  evaluator<typename Derived::Coefficients> ae(_a.coeffs());
94  evaluator<typename OtherDerived::Coefficients> be(_b.coeffs());
95 
96  const double* a = _a.coeffs().data();
97  Packet2d b_xy = be.template packet<BAlignment, Packet2d>(0);
98  Packet2d b_zw = be.template packet<BAlignment, Packet2d>(2);
99  Packet2d a_xx = pset1<Packet2d>(a[0]);
100  Packet2d a_yy = pset1<Packet2d>(a[1]);
101  Packet2d a_zz = pset1<Packet2d>(a[2]);
102  Packet2d a_ww = pset1<Packet2d>(a[3]);
103 
104  // two temporaries:
105  Packet2d t1, t2;
106 
107  /*
108  * t1 = ww*xy + yy*zw
109  * t2 = zz*xy - xx*zw
110  * res.xy = t1 +/- swap(t2)
111  */
112  t1 = padd(pmul(a_ww, b_xy), pmul(a_yy, b_zw));
113  t2 = psub(pmul(a_zz, b_xy), pmul(a_xx, b_zw));
114  pstoret<double, Packet2d, ResAlignment>(&res.x(), paddsub(t1, preverse(t2)));
115 
116  /*
117  * t1 = ww*zw - yy*xy
118  * t2 = zz*zw + xx*xy
119  * res.zw = t1 -/+ swap(t2) = swap( swap(t1) +/- t2)
120  */
121  t1 = psub(pmul(a_ww, b_zw), pmul(a_yy, b_xy));
122  t2 = padd(pmul(a_zz, b_zw), pmul(a_xx, b_xy));
123  pstoret<double, Packet2d, ResAlignment>(&res.z(), preverse(paddsub(preverse(t1), t2)));
124 
125  return res;
126  }
127 };
128 
129 template <class Derived>
130 struct quat_conj<Architecture::Target, Derived, double> {
131  enum { ResAlignment = traits<Quaternion<double> >::Alignment };
132  static inline Quaternion<double> run(const QuaternionBase<Derived>& q) {
133  evaluator<typename Derived::Coefficients> qe(q.coeffs());
134  Quaternion<double> res;
135  const double neg_zero = numext::bit_cast<double>(0x8000000000000000ull);
136  const double arr1[2] = {neg_zero, neg_zero};
137  const double arr2[2] = {neg_zero, 0.0};
138  const Packet2d mask0 = ploadu<Packet2d>(arr1);
139  const Packet2d mask2 = ploadu<Packet2d>(arr2);
140  pstoret<double, Packet2d, ResAlignment>(&res.x(),
141  pxor(mask0, qe.template packet<traits<Derived>::Alignment, Packet2d>(0)));
142  pstoret<double, Packet2d, ResAlignment>(&res.z(),
143  pxor(mask2, qe.template packet<traits<Derived>::Alignment, Packet2d>(2)));
144  return res;
145  }
146 };
147 
148 #endif // end EIGEN_VECTORIZE_SSE_OR_EIGEN_ARCH_ARM64
149 
150 } // end namespace internal
151 
152 } // end namespace Eigen
153 
154 #endif // EIGEN_GEOMETRY_SIMD_H
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1