1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
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 #include "main.h"
12
13 // using namespace Eigen;
14
15 namespace Eigen {
16 namespace internal {
negate(const T & x)17 template<typename T> T negate(const T& x) { return -x; }
18 }
19 }
20
isApproxAbs(const Scalar & a,const Scalar & b,const typename NumTraits<Scalar>::Real & refvalue)21 template<typename Scalar> bool isApproxAbs(const Scalar& a, const Scalar& b, const typename NumTraits<Scalar>::Real& refvalue)
22 {
23 return internal::isMuchSmallerThan(a-b, refvalue);
24 }
25
areApproxAbs(const Scalar * a,const Scalar * b,int size,const typename NumTraits<Scalar>::Real & refvalue)26 template<typename Scalar> bool areApproxAbs(const Scalar* a, const Scalar* b, int size, const typename NumTraits<Scalar>::Real& refvalue)
27 {
28 for (int i=0; i<size; ++i)
29 {
30 if (!isApproxAbs(a[i],b[i],refvalue))
31 {
32 std::cout << "[" << Map<const Matrix<Scalar,1,Dynamic> >(a,size) << "]" << " != " << Map<const Matrix<Scalar,1,Dynamic> >(b,size) << "\n";
33 return false;
34 }
35 }
36 return true;
37 }
38
areApprox(const Scalar * a,const Scalar * b,int size)39 template<typename Scalar> bool areApprox(const Scalar* a, const Scalar* b, int size)
40 {
41 for (int i=0; i<size; ++i)
42 {
43 if (a[i]!=b[i] && !internal::isApprox(a[i],b[i]))
44 {
45 std::cout << "[" << Map<const Matrix<Scalar,1,Dynamic> >(a,size) << "]" << " != " << Map<const Matrix<Scalar,1,Dynamic> >(b,size) << "\n";
46 return false;
47 }
48 }
49 return true;
50 }
51
52
53 #define CHECK_CWISE2(REFOP, POP) { \
54 for (int i=0; i<PacketSize; ++i) \
55 ref[i] = REFOP(data1[i], data1[i+PacketSize]); \
56 internal::pstore(data2, POP(internal::pload<Packet>(data1), internal::pload<Packet>(data1+PacketSize))); \
57 VERIFY(areApprox(ref, data2, PacketSize) && #POP); \
58 }
59
60 #define CHECK_CWISE1(REFOP, POP) { \
61 for (int i=0; i<PacketSize; ++i) \
62 ref[i] = REFOP(data1[i]); \
63 internal::pstore(data2, POP(internal::pload<Packet>(data1))); \
64 VERIFY(areApprox(ref, data2, PacketSize) && #POP); \
65 }
66
67 template<bool Cond,typename Packet>
68 struct packet_helper
69 {
70 template<typename T>
loadpacket_helper71 inline Packet load(const T* from) const { return internal::pload<Packet>(from); }
72
73 template<typename T>
storepacket_helper74 inline void store(T* to, const Packet& x) const { internal::pstore(to,x); }
75 };
76
77 template<typename Packet>
78 struct packet_helper<false,Packet>
79 {
80 template<typename T>
loadpacket_helper81 inline T load(const T* from) const { return *from; }
82
83 template<typename T>
storepacket_helper84 inline void store(T* to, const T& x) const { *to = x; }
85 };
86
87 #define CHECK_CWISE1_IF(COND, REFOP, POP) if(COND) { \
88 packet_helper<COND,Packet> h; \
89 for (int i=0; i<PacketSize; ++i) \
90 ref[i] = REFOP(data1[i]); \
91 h.store(data2, POP(h.load(data1))); \
92 VERIFY(areApprox(ref, data2, PacketSize) && #POP); \
93 }
94
95 #define REF_ADD(a,b) ((a)+(b))
96 #define REF_SUB(a,b) ((a)-(b))
97 #define REF_MUL(a,b) ((a)*(b))
98 #define REF_DIV(a,b) ((a)/(b))
99
packetmath()100 template<typename Scalar> void packetmath()
101 {
102 using std::abs;
103 typedef typename internal::packet_traits<Scalar>::type Packet;
104 const int PacketSize = internal::packet_traits<Scalar>::size;
105 typedef typename NumTraits<Scalar>::Real RealScalar;
106
107 const int size = PacketSize*4;
108 EIGEN_ALIGN16 Scalar data1[internal::packet_traits<Scalar>::size*4];
109 EIGEN_ALIGN16 Scalar data2[internal::packet_traits<Scalar>::size*4];
110 EIGEN_ALIGN16 Packet packets[PacketSize*2];
111 EIGEN_ALIGN16 Scalar ref[internal::packet_traits<Scalar>::size*4];
112 RealScalar refvalue = 0;
113 for (int i=0; i<size; ++i)
114 {
115 data1[i] = internal::random<Scalar>()/RealScalar(PacketSize);
116 data2[i] = internal::random<Scalar>()/RealScalar(PacketSize);
117 refvalue = (std::max)(refvalue,abs(data1[i]));
118 }
119
120 internal::pstore(data2, internal::pload<Packet>(data1));
121 VERIFY(areApprox(data1, data2, PacketSize) && "aligned load/store");
122
123 for (int offset=0; offset<PacketSize; ++offset)
124 {
125 internal::pstore(data2, internal::ploadu<Packet>(data1+offset));
126 VERIFY(areApprox(data1+offset, data2, PacketSize) && "internal::ploadu");
127 }
128
129 for (int offset=0; offset<PacketSize; ++offset)
130 {
131 internal::pstoreu(data2+offset, internal::pload<Packet>(data1));
132 VERIFY(areApprox(data1, data2+offset, PacketSize) && "internal::pstoreu");
133 }
134
135 for (int offset=0; offset<PacketSize; ++offset)
136 {
137 packets[0] = internal::pload<Packet>(data1);
138 packets[1] = internal::pload<Packet>(data1+PacketSize);
139 if (offset==0) internal::palign<0>(packets[0], packets[1]);
140 else if (offset==1) internal::palign<1>(packets[0], packets[1]);
141 else if (offset==2) internal::palign<2>(packets[0], packets[1]);
142 else if (offset==3) internal::palign<3>(packets[0], packets[1]);
143 internal::pstore(data2, packets[0]);
144
145 for (int i=0; i<PacketSize; ++i)
146 ref[i] = data1[i+offset];
147
148 VERIFY(areApprox(ref, data2, PacketSize) && "internal::palign");
149 }
150
151 CHECK_CWISE2(REF_ADD, internal::padd);
152 CHECK_CWISE2(REF_SUB, internal::psub);
153 CHECK_CWISE2(REF_MUL, internal::pmul);
154 #ifndef EIGEN_VECTORIZE_ALTIVEC
155 if (!internal::is_same<Scalar,int>::value)
156 CHECK_CWISE2(REF_DIV, internal::pdiv);
157 #endif
158 CHECK_CWISE1(internal::negate, internal::pnegate);
159 CHECK_CWISE1(numext::conj, internal::pconj);
160
161 for(int offset=0;offset<3;++offset)
162 {
163 for (int i=0; i<PacketSize; ++i)
164 ref[i] = data1[offset];
165 internal::pstore(data2, internal::pset1<Packet>(data1[offset]));
166 VERIFY(areApprox(ref, data2, PacketSize) && "internal::pset1");
167 }
168
169 VERIFY(internal::isApprox(data1[0], internal::pfirst(internal::pload<Packet>(data1))) && "internal::pfirst");
170
171 if(PacketSize>1)
172 {
173 for(int offset=0;offset<4;++offset)
174 {
175 for(int i=0;i<PacketSize/2;++i)
176 ref[2*i+0] = ref[2*i+1] = data1[offset+i];
177 internal::pstore(data2,internal::ploaddup<Packet>(data1+offset));
178 VERIFY(areApprox(ref, data2, PacketSize) && "ploaddup");
179 }
180 }
181
182 ref[0] = 0;
183 for (int i=0; i<PacketSize; ++i)
184 ref[0] += data1[i];
185 VERIFY(isApproxAbs(ref[0], internal::predux(internal::pload<Packet>(data1)), refvalue) && "internal::predux");
186
187 ref[0] = 1;
188 for (int i=0; i<PacketSize; ++i)
189 ref[0] *= data1[i];
190 VERIFY(internal::isApprox(ref[0], internal::predux_mul(internal::pload<Packet>(data1))) && "internal::predux_mul");
191
192 for (int j=0; j<PacketSize; ++j)
193 {
194 ref[j] = 0;
195 for (int i=0; i<PacketSize; ++i)
196 ref[j] += data1[i+j*PacketSize];
197 packets[j] = internal::pload<Packet>(data1+j*PacketSize);
198 }
199 internal::pstore(data2, internal::preduxp(packets));
200 VERIFY(areApproxAbs(ref, data2, PacketSize, refvalue) && "internal::preduxp");
201
202 for (int i=0; i<PacketSize; ++i)
203 ref[i] = data1[PacketSize-i-1];
204 internal::pstore(data2, internal::preverse(internal::pload<Packet>(data1)));
205 VERIFY(areApprox(ref, data2, PacketSize) && "internal::preverse");
206 }
207
packetmath_real()208 template<typename Scalar> void packetmath_real()
209 {
210 using std::abs;
211 typedef typename internal::packet_traits<Scalar>::type Packet;
212 const int PacketSize = internal::packet_traits<Scalar>::size;
213
214 const int size = PacketSize*4;
215 EIGEN_ALIGN16 Scalar data1[internal::packet_traits<Scalar>::size*4];
216 EIGEN_ALIGN16 Scalar data2[internal::packet_traits<Scalar>::size*4];
217 EIGEN_ALIGN16 Scalar ref[internal::packet_traits<Scalar>::size*4];
218
219 for (int i=0; i<size; ++i)
220 {
221 data1[i] = internal::random<Scalar>(-1,1) * std::pow(Scalar(10), internal::random<Scalar>(-3,3));
222 data2[i] = internal::random<Scalar>(-1,1) * std::pow(Scalar(10), internal::random<Scalar>(-3,3));
223 }
224 CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasSin, std::sin, internal::psin);
225 CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasCos, std::cos, internal::pcos);
226 CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasTan, std::tan, internal::ptan);
227
228 for (int i=0; i<size; ++i)
229 {
230 data1[i] = internal::random<Scalar>(-1,1);
231 data2[i] = internal::random<Scalar>(-1,1);
232 }
233 CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasASin, std::asin, internal::pasin);
234 CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasACos, std::acos, internal::pacos);
235
236 for (int i=0; i<size; ++i)
237 {
238 data1[i] = internal::random<Scalar>(-87,88);
239 data2[i] = internal::random<Scalar>(-87,88);
240 }
241 CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasExp, std::exp, internal::pexp);
242
243 for (int i=0; i<size; ++i)
244 {
245 data1[i] = internal::random<Scalar>(0,1) * std::pow(Scalar(10), internal::random<Scalar>(-6,6));
246 data2[i] = internal::random<Scalar>(0,1) * std::pow(Scalar(10), internal::random<Scalar>(-6,6));
247 }
248 if(internal::random<float>(0,1)<0.1)
249 data1[internal::random<int>(0, PacketSize)] = 0;
250 CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasLog, std::log, internal::plog);
251 CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasSqrt, std::sqrt, internal::psqrt);
252 }
253
packetmath_notcomplex()254 template<typename Scalar> void packetmath_notcomplex()
255 {
256 using std::abs;
257 typedef typename internal::packet_traits<Scalar>::type Packet;
258 const int PacketSize = internal::packet_traits<Scalar>::size;
259
260 EIGEN_ALIGN16 Scalar data1[internal::packet_traits<Scalar>::size*4];
261 EIGEN_ALIGN16 Scalar data2[internal::packet_traits<Scalar>::size*4];
262 EIGEN_ALIGN16 Scalar ref[internal::packet_traits<Scalar>::size*4];
263
264 Array<Scalar,Dynamic,1>::Map(data1, internal::packet_traits<Scalar>::size*4).setRandom();
265
266 ref[0] = data1[0];
267 for (int i=0; i<PacketSize; ++i)
268 ref[0] = (std::min)(ref[0],data1[i]);
269 VERIFY(internal::isApprox(ref[0], internal::predux_min(internal::pload<Packet>(data1))) && "internal::predux_min");
270
271 CHECK_CWISE2((std::min), internal::pmin);
272 CHECK_CWISE2((std::max), internal::pmax);
273 CHECK_CWISE1(abs, internal::pabs);
274
275 ref[0] = data1[0];
276 for (int i=0; i<PacketSize; ++i)
277 ref[0] = (std::max)(ref[0],data1[i]);
278 VERIFY(internal::isApprox(ref[0], internal::predux_max(internal::pload<Packet>(data1))) && "internal::predux_max");
279
280 for (int i=0; i<PacketSize; ++i)
281 ref[i] = data1[0]+Scalar(i);
282 internal::pstore(data2, internal::plset(data1[0]));
283 VERIFY(areApprox(ref, data2, PacketSize) && "internal::plset");
284 }
285
test_conj_helper(Scalar * data1,Scalar * data2,Scalar * ref,Scalar * pval)286 template<typename Scalar,bool ConjLhs,bool ConjRhs> void test_conj_helper(Scalar* data1, Scalar* data2, Scalar* ref, Scalar* pval)
287 {
288 typedef typename internal::packet_traits<Scalar>::type Packet;
289 const int PacketSize = internal::packet_traits<Scalar>::size;
290
291 internal::conj_if<ConjLhs> cj0;
292 internal::conj_if<ConjRhs> cj1;
293 internal::conj_helper<Scalar,Scalar,ConjLhs,ConjRhs> cj;
294 internal::conj_helper<Packet,Packet,ConjLhs,ConjRhs> pcj;
295
296 for(int i=0;i<PacketSize;++i)
297 {
298 ref[i] = cj0(data1[i]) * cj1(data2[i]);
299 VERIFY(internal::isApprox(ref[i], cj.pmul(data1[i],data2[i])) && "conj_helper pmul");
300 }
301 internal::pstore(pval,pcj.pmul(internal::pload<Packet>(data1),internal::pload<Packet>(data2)));
302 VERIFY(areApprox(ref, pval, PacketSize) && "conj_helper pmul");
303
304 for(int i=0;i<PacketSize;++i)
305 {
306 Scalar tmp = ref[i];
307 ref[i] += cj0(data1[i]) * cj1(data2[i]);
308 VERIFY(internal::isApprox(ref[i], cj.pmadd(data1[i],data2[i],tmp)) && "conj_helper pmadd");
309 }
310 internal::pstore(pval,pcj.pmadd(internal::pload<Packet>(data1),internal::pload<Packet>(data2),internal::pload<Packet>(pval)));
311 VERIFY(areApprox(ref, pval, PacketSize) && "conj_helper pmadd");
312 }
313
packetmath_complex()314 template<typename Scalar> void packetmath_complex()
315 {
316 typedef typename internal::packet_traits<Scalar>::type Packet;
317 const int PacketSize = internal::packet_traits<Scalar>::size;
318
319 const int size = PacketSize*4;
320 EIGEN_ALIGN16 Scalar data1[PacketSize*4];
321 EIGEN_ALIGN16 Scalar data2[PacketSize*4];
322 EIGEN_ALIGN16 Scalar ref[PacketSize*4];
323 EIGEN_ALIGN16 Scalar pval[PacketSize*4];
324
325 for (int i=0; i<size; ++i)
326 {
327 data1[i] = internal::random<Scalar>() * Scalar(1e2);
328 data2[i] = internal::random<Scalar>() * Scalar(1e2);
329 }
330
331 test_conj_helper<Scalar,false,false> (data1,data2,ref,pval);
332 test_conj_helper<Scalar,false,true> (data1,data2,ref,pval);
333 test_conj_helper<Scalar,true,false> (data1,data2,ref,pval);
334 test_conj_helper<Scalar,true,true> (data1,data2,ref,pval);
335
336 {
337 for(int i=0;i<PacketSize;++i)
338 ref[i] = Scalar(std::imag(data1[i]),std::real(data1[i]));
339 internal::pstore(pval,internal::pcplxflip(internal::pload<Packet>(data1)));
340 VERIFY(areApprox(ref, pval, PacketSize) && "pcplxflip");
341 }
342
343
344 }
345
test_packetmath()346 void test_packetmath()
347 {
348 for(int i = 0; i < g_repeat; i++) {
349 CALL_SUBTEST_1( packetmath<float>() );
350 CALL_SUBTEST_2( packetmath<double>() );
351 CALL_SUBTEST_3( packetmath<int>() );
352 CALL_SUBTEST_1( packetmath<std::complex<float> >() );
353 CALL_SUBTEST_2( packetmath<std::complex<double> >() );
354
355 CALL_SUBTEST_1( packetmath_notcomplex<float>() );
356 CALL_SUBTEST_2( packetmath_notcomplex<double>() );
357 CALL_SUBTEST_3( packetmath_notcomplex<int>() );
358
359 CALL_SUBTEST_1( packetmath_real<float>() );
360 CALL_SUBTEST_2( packetmath_real<double>() );
361
362 CALL_SUBTEST_1( packetmath_complex<std::complex<float> >() );
363 CALL_SUBTEST_2( packetmath_complex<std::complex<double> >() );
364 }
365 }
366