Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-22 10:34:17

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2009 Rohit Garg <rpg.314@gmail.com>
0005 // Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
0006 //
0007 // This Source Code Form is subject to the terms of the Mozilla
0008 // Public License v. 2.0. If a copy of the MPL was not distributed
0009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0010 
0011 #ifndef EIGEN_GEOMETRY_SIMD_H
0012 #define EIGEN_GEOMETRY_SIMD_H
0013 
0014 namespace Eigen { 
0015 
0016 namespace internal {
0017 
0018 template<class Derived, class OtherDerived>
0019 struct quat_product<Architecture::Target, Derived, OtherDerived, float>
0020 {
0021   enum {
0022     AAlignment = traits<Derived>::Alignment,
0023     BAlignment = traits<OtherDerived>::Alignment,
0024     ResAlignment = traits<Quaternion<float> >::Alignment
0025   };
0026   static inline Quaternion<float> run(const QuaternionBase<Derived>& _a, const QuaternionBase<OtherDerived>& _b)
0027   {
0028     evaluator<typename Derived::Coefficients> ae(_a.coeffs());
0029     evaluator<typename OtherDerived::Coefficients> be(_b.coeffs());
0030     Quaternion<float> res;
0031     const float neg_zero = numext::bit_cast<float>(0x80000000u);
0032     const float arr[4] = {0.f, 0.f, 0.f, neg_zero};
0033     const Packet4f mask = ploadu<Packet4f>(arr);
0034     Packet4f a = ae.template packet<AAlignment,Packet4f>(0);
0035     Packet4f b = be.template packet<BAlignment,Packet4f>(0);
0036     Packet4f s1 = pmul(vec4f_swizzle1(a,1,2,0,2),vec4f_swizzle1(b,2,0,1,2));
0037     Packet4f s2 = pmul(vec4f_swizzle1(a,3,3,3,1),vec4f_swizzle1(b,0,1,2,1));
0038     pstoret<float,Packet4f,ResAlignment>(
0039               &res.x(),
0040               padd(psub(pmul(a,vec4f_swizzle1(b,3,3,3,3)),
0041                                     pmul(vec4f_swizzle1(a,2,0,1,0),
0042                                                vec4f_swizzle1(b,1,2,0,0))),
0043                          pxor(mask,padd(s1,s2))));
0044     
0045     return res;
0046   }
0047 };
0048 
0049 template<class Derived>
0050 struct quat_conj<Architecture::Target, Derived, float>
0051 {
0052   enum {
0053     ResAlignment = traits<Quaternion<float> >::Alignment
0054   };
0055   static inline Quaternion<float> run(const QuaternionBase<Derived>& q)
0056   {
0057     evaluator<typename Derived::Coefficients> qe(q.coeffs());
0058     Quaternion<float> res;
0059     const float neg_zero = numext::bit_cast<float>(0x80000000u);
0060     const float arr[4] = {neg_zero, neg_zero, neg_zero,0.f};
0061     const Packet4f mask = ploadu<Packet4f>(arr);
0062     pstoret<float,Packet4f,ResAlignment>(&res.x(), pxor(mask, qe.template packet<traits<Derived>::Alignment,Packet4f>(0)));
0063     return res;
0064   }
0065 };
0066 
0067 
0068 template<typename VectorLhs,typename VectorRhs>
0069 struct cross3_impl<Architecture::Target,VectorLhs,VectorRhs,float,true>
0070 {
0071   enum {
0072     ResAlignment = traits<typename plain_matrix_type<VectorLhs>::type>::Alignment
0073   };
0074   static inline typename plain_matrix_type<VectorLhs>::type
0075   run(const VectorLhs& lhs, const VectorRhs& rhs)
0076   {
0077     evaluator<VectorLhs> lhs_eval(lhs);
0078     evaluator<VectorRhs> rhs_eval(rhs);
0079     Packet4f a = lhs_eval.template packet<traits<VectorLhs>::Alignment,Packet4f>(0);
0080     Packet4f b = rhs_eval.template packet<traits<VectorRhs>::Alignment,Packet4f>(0);
0081     Packet4f mul1 = pmul(vec4f_swizzle1(a,1,2,0,3),vec4f_swizzle1(b,2,0,1,3));
0082     Packet4f mul2 = pmul(vec4f_swizzle1(a,2,0,1,3),vec4f_swizzle1(b,1,2,0,3));
0083     typename plain_matrix_type<VectorLhs>::type res;
0084     pstoret<float,Packet4f,ResAlignment>(&res.x(),psub(mul1,mul2));
0085     return res;
0086   }
0087 };
0088 
0089 
0090 
0091 #if (defined EIGEN_VECTORIZE_SSE) || (EIGEN_ARCH_ARM64)
0092 
0093 template<class Derived, class OtherDerived>
0094 struct quat_product<Architecture::Target, Derived, OtherDerived, double>
0095 {
0096   enum {
0097     BAlignment = traits<OtherDerived>::Alignment,
0098     ResAlignment = traits<Quaternion<double> >::Alignment
0099   };
0100 
0101   static inline Quaternion<double> run(const QuaternionBase<Derived>& _a, const QuaternionBase<OtherDerived>& _b)
0102   {
0103   Quaternion<double> res;
0104 
0105   evaluator<typename Derived::Coefficients> ae(_a.coeffs());
0106   evaluator<typename OtherDerived::Coefficients> be(_b.coeffs());
0107 
0108   const double* a = _a.coeffs().data();
0109   Packet2d b_xy = be.template packet<BAlignment,Packet2d>(0);
0110   Packet2d b_zw = be.template packet<BAlignment,Packet2d>(2);
0111   Packet2d a_xx = pset1<Packet2d>(a[0]);
0112   Packet2d a_yy = pset1<Packet2d>(a[1]);
0113   Packet2d a_zz = pset1<Packet2d>(a[2]);
0114   Packet2d a_ww = pset1<Packet2d>(a[3]);
0115 
0116   // two temporaries:
0117   Packet2d t1, t2;
0118 
0119   /*
0120    * t1 = ww*xy + yy*zw
0121    * t2 = zz*xy - xx*zw
0122    * res.xy = t1 +/- swap(t2)
0123    */
0124   t1 = padd(pmul(a_ww, b_xy), pmul(a_yy, b_zw));
0125   t2 = psub(pmul(a_zz, b_xy), pmul(a_xx, b_zw));
0126   pstoret<double,Packet2d,ResAlignment>(&res.x(), paddsub(t1, preverse(t2)));
0127   
0128   /*
0129    * t1 = ww*zw - yy*xy
0130    * t2 = zz*zw + xx*xy
0131    * res.zw = t1 -/+ swap(t2) = swap( swap(t1) +/- t2)
0132    */
0133   t1 = psub(pmul(a_ww, b_zw), pmul(a_yy, b_xy));
0134   t2 = padd(pmul(a_zz, b_zw), pmul(a_xx, b_xy));
0135   pstoret<double,Packet2d,ResAlignment>(&res.z(), preverse(paddsub(preverse(t1), t2)));
0136 
0137   return res;
0138 }
0139 };
0140 
0141 template<class Derived>
0142 struct quat_conj<Architecture::Target, Derived, double>
0143 {
0144   enum {
0145     ResAlignment = traits<Quaternion<double> >::Alignment
0146   };
0147   static inline Quaternion<double> run(const QuaternionBase<Derived>& q)
0148   {
0149     evaluator<typename Derived::Coefficients> qe(q.coeffs());
0150     Quaternion<double> res;
0151     const double neg_zero = numext::bit_cast<double>(0x8000000000000000ull);
0152     const double arr1[2] = {neg_zero, neg_zero};
0153     const double arr2[2] = {neg_zero,  0.0};
0154     const Packet2d mask0 = ploadu<Packet2d>(arr1);
0155     const Packet2d mask2 = ploadu<Packet2d>(arr2);
0156     pstoret<double,Packet2d,ResAlignment>(&res.x(), pxor(mask0, qe.template packet<traits<Derived>::Alignment,Packet2d>(0)));
0157     pstoret<double,Packet2d,ResAlignment>(&res.z(), pxor(mask2, qe.template packet<traits<Derived>::Alignment,Packet2d>(2)));
0158     return res;
0159   }
0160 };
0161 
0162 #endif // end EIGEN_VECTORIZE_SSE_OR_EIGEN_ARCH_ARM64
0163 
0164 } // end namespace internal
0165 
0166 } // end namespace Eigen
0167 
0168 #endif // EIGEN_GEOMETRY_SIMD_H