Complex.h revision 7faaa9f3f0df9d23790277834d426c3d992ac3ba
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2010 Gael Guennebaud <gael.guennebaud@inria.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_COMPLEX_NEON_H
11#define EIGEN_COMPLEX_NEON_H
12
13namespace Eigen {
14
15namespace internal {
16
17static uint32x4_t p4ui_CONJ_XOR = EIGEN_INIT_NEON_PACKET4(0x00000000, 0x80000000, 0x00000000, 0x80000000);
18static uint32x2_t p2ui_CONJ_XOR = EIGEN_INIT_NEON_PACKET2(0x00000000, 0x80000000);
19
20//---------- float ----------
21struct Packet2cf
22{
23  EIGEN_STRONG_INLINE Packet2cf() {}
24  EIGEN_STRONG_INLINE explicit Packet2cf(const Packet4f& a) : v(a) {}
25  Packet4f  v;
26};
27
28template<> struct packet_traits<std::complex<float> >  : default_packet_traits
29{
30  typedef Packet2cf type;
31  enum {
32    Vectorizable = 1,
33    AlignedOnScalar = 1,
34    size = 2,
35
36    HasAdd    = 1,
37    HasSub    = 1,
38    HasMul    = 1,
39    HasDiv    = 1,
40    HasNegate = 1,
41    HasAbs    = 0,
42    HasAbs2   = 0,
43    HasMin    = 0,
44    HasMax    = 0,
45    HasSetLinear = 0
46  };
47};
48
49template<> struct unpacket_traits<Packet2cf> { typedef std::complex<float> type; enum {size=2}; };
50
51template<> EIGEN_STRONG_INLINE Packet2cf pset1<Packet2cf>(const std::complex<float>&  from)
52{
53  float32x2_t r64;
54  r64 = vld1_f32((float *)&from);
55
56  return Packet2cf(vcombine_f32(r64, r64));
57}
58
59template<> EIGEN_STRONG_INLINE Packet2cf padd<Packet2cf>(const Packet2cf& a, const Packet2cf& b) { return Packet2cf(padd<Packet4f>(a.v,b.v)); }
60template<> EIGEN_STRONG_INLINE Packet2cf psub<Packet2cf>(const Packet2cf& a, const Packet2cf& b) { return Packet2cf(psub<Packet4f>(a.v,b.v)); }
61template<> EIGEN_STRONG_INLINE Packet2cf pnegate(const Packet2cf& a) { return Packet2cf(pnegate<Packet4f>(a.v)); }
62template<> EIGEN_STRONG_INLINE Packet2cf pconj(const Packet2cf& a)
63{
64  Packet4ui b = vreinterpretq_u32_f32(a.v);
65  return Packet2cf(vreinterpretq_f32_u32(veorq_u32(b, p4ui_CONJ_XOR)));
66}
67
68template<> EIGEN_STRONG_INLINE Packet2cf pmul<Packet2cf>(const Packet2cf& a, const Packet2cf& b)
69{
70  Packet4f v1, v2;
71
72  // Get the real values of a | a1_re | a1_re | a2_re | a2_re |
73  v1 = vcombine_f32(vdup_lane_f32(vget_low_f32(a.v), 0), vdup_lane_f32(vget_high_f32(a.v), 0));
74  // Get the real values of a | a1_im | a1_im | a2_im | a2_im |
75  v2 = vcombine_f32(vdup_lane_f32(vget_low_f32(a.v), 1), vdup_lane_f32(vget_high_f32(a.v), 1));
76  // Multiply the real a with b
77  v1 = vmulq_f32(v1, b.v);
78  // Multiply the imag a with b
79  v2 = vmulq_f32(v2, b.v);
80  // Conjugate v2
81  v2 = vreinterpretq_f32_u32(veorq_u32(vreinterpretq_u32_f32(v2), p4ui_CONJ_XOR));
82  // Swap real/imag elements in v2.
83  v2 = vrev64q_f32(v2);
84  // Add and return the result
85  return Packet2cf(vaddq_f32(v1, v2));
86}
87
88template<> EIGEN_STRONG_INLINE Packet2cf pand   <Packet2cf>(const Packet2cf& a, const Packet2cf& b)
89{
90  return Packet2cf(vreinterpretq_f32_u32(vorrq_u32(vreinterpretq_u32_f32(a.v),vreinterpretq_u32_f32(b.v))));
91}
92template<> EIGEN_STRONG_INLINE Packet2cf por    <Packet2cf>(const Packet2cf& a, const Packet2cf& b)
93{
94  return Packet2cf(vreinterpretq_f32_u32(vorrq_u32(vreinterpretq_u32_f32(a.v),vreinterpretq_u32_f32(b.v))));
95}
96template<> EIGEN_STRONG_INLINE Packet2cf pxor   <Packet2cf>(const Packet2cf& a, const Packet2cf& b)
97{
98  return Packet2cf(vreinterpretq_f32_u32(veorq_u32(vreinterpretq_u32_f32(a.v),vreinterpretq_u32_f32(b.v))));
99}
100template<> EIGEN_STRONG_INLINE Packet2cf pandnot<Packet2cf>(const Packet2cf& a, const Packet2cf& b)
101{
102  return Packet2cf(vreinterpretq_f32_u32(vbicq_u32(vreinterpretq_u32_f32(a.v),vreinterpretq_u32_f32(b.v))));
103}
104
105template<> EIGEN_STRONG_INLINE Packet2cf pload<Packet2cf>(const std::complex<float>* from) { EIGEN_DEBUG_ALIGNED_LOAD return Packet2cf(pload<Packet4f>((const float*)from)); }
106template<> EIGEN_STRONG_INLINE Packet2cf ploadu<Packet2cf>(const std::complex<float>* from) { EIGEN_DEBUG_UNALIGNED_LOAD return Packet2cf(ploadu<Packet4f>((const float*)from)); }
107
108template<> EIGEN_STRONG_INLINE Packet2cf ploaddup<Packet2cf>(const std::complex<float>* from) { return pset1<Packet2cf>(*from); }
109
110template<> EIGEN_STRONG_INLINE void pstore <std::complex<float> >(std::complex<float> *   to, const Packet2cf& from) { EIGEN_DEBUG_ALIGNED_STORE pstore((float*)to, from.v); }
111template<> EIGEN_STRONG_INLINE void pstoreu<std::complex<float> >(std::complex<float> *   to, const Packet2cf& from) { EIGEN_DEBUG_UNALIGNED_STORE pstoreu((float*)to, from.v); }
112
113template<> EIGEN_STRONG_INLINE void prefetch<std::complex<float> >(const std::complex<float> *   addr) { __pld((float *)addr); }
114
115template<> EIGEN_STRONG_INLINE std::complex<float>  pfirst<Packet2cf>(const Packet2cf& a)
116{
117  std::complex<float> EIGEN_ALIGN16 x[2];
118  vst1q_f32((float *)x, a.v);
119  return x[0];
120}
121
122template<> EIGEN_STRONG_INLINE Packet2cf preverse(const Packet2cf& a)
123{
124  float32x2_t a_lo, a_hi;
125  Packet4f a_r128;
126
127  a_lo = vget_low_f32(a.v);
128  a_hi = vget_high_f32(a.v);
129  a_r128 = vcombine_f32(a_hi, a_lo);
130
131  return Packet2cf(a_r128);
132}
133
134template<> EIGEN_STRONG_INLINE Packet2cf pcplxflip<Packet2cf>(const Packet2cf& a)
135{
136  return Packet2cf(vrev64q_f32(a.v));
137}
138
139template<> EIGEN_STRONG_INLINE std::complex<float> predux<Packet2cf>(const Packet2cf& a)
140{
141  float32x2_t a1, a2;
142  std::complex<float> s;
143
144  a1 = vget_low_f32(a.v);
145  a2 = vget_high_f32(a.v);
146  a2 = vadd_f32(a1, a2);
147  vst1_f32((float *)&s, a2);
148
149  return s;
150}
151
152template<> EIGEN_STRONG_INLINE Packet2cf preduxp<Packet2cf>(const Packet2cf* vecs)
153{
154  Packet4f sum1, sum2, sum;
155
156  // Add the first two 64-bit float32x2_t of vecs[0]
157  sum1 = vcombine_f32(vget_low_f32(vecs[0].v), vget_low_f32(vecs[1].v));
158  sum2 = vcombine_f32(vget_high_f32(vecs[0].v), vget_high_f32(vecs[1].v));
159  sum = vaddq_f32(sum1, sum2);
160
161  return Packet2cf(sum);
162}
163
164template<> EIGEN_STRONG_INLINE std::complex<float> predux_mul<Packet2cf>(const Packet2cf& a)
165{
166  float32x2_t a1, a2, v1, v2, prod;
167  std::complex<float> s;
168
169  a1 = vget_low_f32(a.v);
170  a2 = vget_high_f32(a.v);
171   // Get the real values of a | a1_re | a1_re | a2_re | a2_re |
172  v1 = vdup_lane_f32(a1, 0);
173  // Get the real values of a | a1_im | a1_im | a2_im | a2_im |
174  v2 = vdup_lane_f32(a1, 1);
175  // Multiply the real a with b
176  v1 = vmul_f32(v1, a2);
177  // Multiply the imag a with b
178  v2 = vmul_f32(v2, a2);
179  // Conjugate v2
180  v2 = vreinterpret_f32_u32(veor_u32(vreinterpret_u32_f32(v2), p2ui_CONJ_XOR));
181  // Swap real/imag elements in v2.
182  v2 = vrev64_f32(v2);
183  // Add v1, v2
184  prod = vadd_f32(v1, v2);
185
186  vst1_f32((float *)&s, prod);
187
188  return s;
189}
190
191template<int Offset>
192struct palign_impl<Offset,Packet2cf>
193{
194  EIGEN_STRONG_INLINE static void run(Packet2cf& first, const Packet2cf& second)
195  {
196    if (Offset==1)
197    {
198      first.v = vextq_f32(first.v, second.v, 2);
199    }
200  }
201};
202
203template<> struct conj_helper<Packet2cf, Packet2cf, false,true>
204{
205  EIGEN_STRONG_INLINE Packet2cf pmadd(const Packet2cf& x, const Packet2cf& y, const Packet2cf& c) const
206  { return padd(pmul(x,y),c); }
207
208  EIGEN_STRONG_INLINE Packet2cf pmul(const Packet2cf& a, const Packet2cf& b) const
209  {
210    return internal::pmul(a, pconj(b));
211  }
212};
213
214template<> struct conj_helper<Packet2cf, Packet2cf, true,false>
215{
216  EIGEN_STRONG_INLINE Packet2cf pmadd(const Packet2cf& x, const Packet2cf& y, const Packet2cf& c) const
217  { return padd(pmul(x,y),c); }
218
219  EIGEN_STRONG_INLINE Packet2cf pmul(const Packet2cf& a, const Packet2cf& b) const
220  {
221    return internal::pmul(pconj(a), b);
222  }
223};
224
225template<> struct conj_helper<Packet2cf, Packet2cf, true,true>
226{
227  EIGEN_STRONG_INLINE Packet2cf pmadd(const Packet2cf& x, const Packet2cf& y, const Packet2cf& c) const
228  { return padd(pmul(x,y),c); }
229
230  EIGEN_STRONG_INLINE Packet2cf pmul(const Packet2cf& a, const Packet2cf& b) const
231  {
232    return pconj(internal::pmul(a, b));
233  }
234};
235
236template<> EIGEN_STRONG_INLINE Packet2cf pdiv<Packet2cf>(const Packet2cf& a, const Packet2cf& b)
237{
238  // TODO optimize it for AltiVec
239  Packet2cf res = conj_helper<Packet2cf,Packet2cf,false,true>().pmul(a,b);
240  Packet4f s, rev_s;
241
242  // this computes the norm
243  s = vmulq_f32(b.v, b.v);
244  rev_s = vrev64q_f32(s);
245
246  return Packet2cf(pdiv(res.v, vaddq_f32(s,rev_s)));
247}
248
249} // end namespace internal
250
251} // end namespace Eigen
252
253#endif // EIGEN_COMPLEX_NEON_H
254