2008-08-21 04:08:38 +08:00
|
|
|
// This file is part of Eigen, a lightweight C++ template library
|
2009-05-23 02:25:33 +08:00
|
|
|
// for linear algebra.
|
2008-08-21 04:08:38 +08:00
|
|
|
//
|
2010-06-25 05:21:58 +08:00
|
|
|
// Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
|
2008-11-24 21:40:43 +08:00
|
|
|
// Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
|
2008-08-21 04:08:38 +08:00
|
|
|
//
|
2012-07-14 02:42:47 +08:00
|
|
|
// This Source Code Form is subject to the terms of the Mozilla
|
|
|
|
// Public License v. 2.0. If a copy of the MPL was not distributed
|
|
|
|
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
2008-08-21 04:08:38 +08:00
|
|
|
|
|
|
|
#include "main.h"
|
|
|
|
|
|
|
|
// using namespace Eigen;
|
|
|
|
|
2010-10-25 22:15:22 +08:00
|
|
|
namespace Eigen {
|
|
|
|
namespace internal {
|
|
|
|
template<typename T> T negate(const T& x) { return -x; }
|
|
|
|
}
|
|
|
|
}
|
2009-03-20 18:03:24 +08:00
|
|
|
|
2010-07-05 16:54:24 +08:00
|
|
|
template<typename Scalar> bool isApproxAbs(const Scalar& a, const Scalar& b, const typename NumTraits<Scalar>::Real& refvalue)
|
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
return internal::isMuchSmallerThan(a-b, refvalue);
|
2010-07-05 16:54:24 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
template<typename Scalar> bool areApproxAbs(const Scalar* a, const Scalar* b, int size, const typename NumTraits<Scalar>::Real& refvalue)
|
|
|
|
{
|
|
|
|
for (int i=0; i<size; ++i)
|
|
|
|
{
|
|
|
|
if (!isApproxAbs(a[i],b[i],refvalue))
|
|
|
|
{
|
2011-05-19 03:11:03 +08:00
|
|
|
std::cout << "[" << Map<const Matrix<Scalar,1,Dynamic> >(a,size) << "]" << " != " << Map<const Matrix<Scalar,1,Dynamic> >(b,size) << "\n";
|
2010-07-05 16:54:24 +08:00
|
|
|
return false;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
2008-08-21 04:08:38 +08:00
|
|
|
template<typename Scalar> bool areApprox(const Scalar* a, const Scalar* b, int size)
|
|
|
|
{
|
|
|
|
for (int i=0; i<size; ++i)
|
2010-07-05 16:54:24 +08:00
|
|
|
{
|
2013-02-15 06:34:05 +08:00
|
|
|
if (a[i]!=b[i] && !internal::isApprox(a[i],b[i]))
|
2010-07-05 16:54:24 +08:00
|
|
|
{
|
2011-05-19 03:11:03 +08:00
|
|
|
std::cout << "[" << Map<const Matrix<Scalar,1,Dynamic> >(a,size) << "]" << " != " << Map<const Matrix<Scalar,1,Dynamic> >(b,size) << "\n";
|
2010-07-05 16:54:24 +08:00
|
|
|
return false;
|
2010-03-04 01:25:41 +08:00
|
|
|
}
|
2010-07-05 16:54:24 +08:00
|
|
|
}
|
2008-08-21 04:08:38 +08:00
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
2010-07-05 16:54:24 +08:00
|
|
|
|
2009-03-10 02:40:09 +08:00
|
|
|
#define CHECK_CWISE2(REFOP, POP) { \
|
2008-08-21 04:08:38 +08:00
|
|
|
for (int i=0; i<PacketSize; ++i) \
|
|
|
|
ref[i] = REFOP(data1[i], data1[i+PacketSize]); \
|
2010-10-25 22:15:22 +08:00
|
|
|
internal::pstore(data2, POP(internal::pload<Packet>(data1), internal::pload<Packet>(data1+PacketSize))); \
|
2008-08-21 04:08:38 +08:00
|
|
|
VERIFY(areApprox(ref, data2, PacketSize) && #POP); \
|
|
|
|
}
|
|
|
|
|
2009-03-10 02:40:09 +08:00
|
|
|
#define CHECK_CWISE1(REFOP, POP) { \
|
|
|
|
for (int i=0; i<PacketSize; ++i) \
|
|
|
|
ref[i] = REFOP(data1[i]); \
|
2010-10-25 22:15:22 +08:00
|
|
|
internal::pstore(data2, POP(internal::pload<Packet>(data1))); \
|
2009-03-10 02:40:09 +08:00
|
|
|
VERIFY(areApprox(ref, data2, PacketSize) && #POP); \
|
|
|
|
}
|
|
|
|
|
2009-03-25 20:26:13 +08:00
|
|
|
template<bool Cond,typename Packet>
|
|
|
|
struct packet_helper
|
|
|
|
{
|
|
|
|
template<typename T>
|
2010-10-25 22:15:22 +08:00
|
|
|
inline Packet load(const T* from) const { return internal::pload<Packet>(from); }
|
2010-07-05 16:54:24 +08:00
|
|
|
|
2009-03-25 20:26:13 +08:00
|
|
|
template<typename T>
|
2010-10-25 22:15:22 +08:00
|
|
|
inline void store(T* to, const Packet& x) const { internal::pstore(to,x); }
|
2009-03-25 20:26:13 +08:00
|
|
|
};
|
|
|
|
|
|
|
|
template<typename Packet>
|
|
|
|
struct packet_helper<false,Packet>
|
|
|
|
{
|
|
|
|
template<typename T>
|
|
|
|
inline T load(const T* from) const { return *from; }
|
2010-07-05 16:54:24 +08:00
|
|
|
|
2009-03-25 20:26:13 +08:00
|
|
|
template<typename T>
|
|
|
|
inline void store(T* to, const T& x) const { *to = x; }
|
|
|
|
};
|
|
|
|
|
|
|
|
#define CHECK_CWISE1_IF(COND, REFOP, POP) if(COND) { \
|
|
|
|
packet_helper<COND,Packet> h; \
|
|
|
|
for (int i=0; i<PacketSize; ++i) \
|
|
|
|
ref[i] = REFOP(data1[i]); \
|
|
|
|
h.store(data2, POP(h.load(data1))); \
|
|
|
|
VERIFY(areApprox(ref, data2, PacketSize) && #POP); \
|
|
|
|
}
|
|
|
|
|
2008-08-21 04:08:38 +08:00
|
|
|
#define REF_ADD(a,b) ((a)+(b))
|
|
|
|
#define REF_SUB(a,b) ((a)-(b))
|
|
|
|
#define REF_MUL(a,b) ((a)*(b))
|
|
|
|
#define REF_DIV(a,b) ((a)/(b))
|
|
|
|
|
|
|
|
template<typename Scalar> void packetmath()
|
|
|
|
{
|
2012-11-06 22:25:50 +08:00
|
|
|
using std::abs;
|
2010-10-25 22:15:22 +08:00
|
|
|
typedef typename internal::packet_traits<Scalar>::type Packet;
|
|
|
|
const int PacketSize = internal::packet_traits<Scalar>::size;
|
2010-07-05 16:54:24 +08:00
|
|
|
typedef typename NumTraits<Scalar>::Real RealScalar;
|
2008-08-21 04:08:38 +08:00
|
|
|
|
|
|
|
const int size = PacketSize*4;
|
2010-10-25 22:15:22 +08:00
|
|
|
EIGEN_ALIGN16 Scalar data1[internal::packet_traits<Scalar>::size*4];
|
|
|
|
EIGEN_ALIGN16 Scalar data2[internal::packet_traits<Scalar>::size*4];
|
2009-10-05 22:11:11 +08:00
|
|
|
EIGEN_ALIGN16 Packet packets[PacketSize*2];
|
2010-10-25 22:15:22 +08:00
|
|
|
EIGEN_ALIGN16 Scalar ref[internal::packet_traits<Scalar>::size*4];
|
2010-07-05 16:54:24 +08:00
|
|
|
RealScalar refvalue = 0;
|
2008-08-21 04:08:38 +08:00
|
|
|
for (int i=0; i<size; ++i)
|
|
|
|
{
|
2011-02-23 23:20:55 +08:00
|
|
|
data1[i] = internal::random<Scalar>()/RealScalar(PacketSize);
|
|
|
|
data2[i] = internal::random<Scalar>()/RealScalar(PacketSize);
|
2012-11-06 22:25:50 +08:00
|
|
|
refvalue = (std::max)(refvalue,abs(data1[i]));
|
2008-08-21 04:08:38 +08:00
|
|
|
}
|
|
|
|
|
2010-10-25 22:15:22 +08:00
|
|
|
internal::pstore(data2, internal::pload<Packet>(data1));
|
2008-08-21 04:08:38 +08:00
|
|
|
VERIFY(areApprox(data1, data2, PacketSize) && "aligned load/store");
|
|
|
|
|
|
|
|
for (int offset=0; offset<PacketSize; ++offset)
|
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
internal::pstore(data2, internal::ploadu<Packet>(data1+offset));
|
|
|
|
VERIFY(areApprox(data1+offset, data2, PacketSize) && "internal::ploadu");
|
2008-08-21 04:08:38 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
for (int offset=0; offset<PacketSize; ++offset)
|
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
internal::pstoreu(data2+offset, internal::pload<Packet>(data1));
|
|
|
|
VERIFY(areApprox(data1, data2+offset, PacketSize) && "internal::pstoreu");
|
2008-08-21 04:08:38 +08:00
|
|
|
}
|
|
|
|
|
2008-08-24 23:27:05 +08:00
|
|
|
for (int offset=0; offset<PacketSize; ++offset)
|
2008-08-21 04:08:38 +08:00
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
packets[0] = internal::pload<Packet>(data1);
|
|
|
|
packets[1] = internal::pload<Packet>(data1+PacketSize);
|
|
|
|
if (offset==0) internal::palign<0>(packets[0], packets[1]);
|
|
|
|
else if (offset==1) internal::palign<1>(packets[0], packets[1]);
|
|
|
|
else if (offset==2) internal::palign<2>(packets[0], packets[1]);
|
|
|
|
else if (offset==3) internal::palign<3>(packets[0], packets[1]);
|
|
|
|
internal::pstore(data2, packets[0]);
|
2008-08-24 23:27:05 +08:00
|
|
|
|
|
|
|
for (int i=0; i<PacketSize; ++i)
|
|
|
|
ref[i] = data1[i+offset];
|
|
|
|
|
2010-10-25 22:15:22 +08:00
|
|
|
VERIFY(areApprox(ref, data2, PacketSize) && "internal::palign");
|
2008-08-21 04:08:38 +08:00
|
|
|
}
|
|
|
|
|
2010-10-25 22:15:22 +08:00
|
|
|
CHECK_CWISE2(REF_ADD, internal::padd);
|
|
|
|
CHECK_CWISE2(REF_SUB, internal::psub);
|
|
|
|
CHECK_CWISE2(REF_MUL, internal::pmul);
|
2008-08-21 22:03:17 +08:00
|
|
|
#ifndef EIGEN_VECTORIZE_ALTIVEC
|
2010-10-26 04:13:49 +08:00
|
|
|
if (!internal::is_same<Scalar,int>::value)
|
2010-10-25 22:15:22 +08:00
|
|
|
CHECK_CWISE2(REF_DIV, internal::pdiv);
|
2008-08-21 22:03:17 +08:00
|
|
|
#endif
|
2010-10-25 22:15:22 +08:00
|
|
|
CHECK_CWISE1(internal::negate, internal::pnegate);
|
|
|
|
CHECK_CWISE1(internal::conj, internal::pconj);
|
2008-08-21 04:08:38 +08:00
|
|
|
|
2011-02-24 05:22:10 +08:00
|
|
|
for(int offset=0;offset<3;++offset)
|
2011-02-24 02:24:26 +08:00
|
|
|
{
|
|
|
|
for (int i=0; i<PacketSize; ++i)
|
|
|
|
ref[i] = data1[offset];
|
|
|
|
internal::pstore(data2, internal::pset1<Packet>(data1[offset]));
|
|
|
|
VERIFY(areApprox(ref, data2, PacketSize) && "internal::pset1");
|
|
|
|
}
|
2011-05-19 03:11:03 +08:00
|
|
|
|
2010-10-25 22:15:22 +08:00
|
|
|
VERIFY(internal::isApprox(data1[0], internal::pfirst(internal::pload<Packet>(data1))) && "internal::pfirst");
|
2011-02-23 23:20:55 +08:00
|
|
|
|
|
|
|
if(PacketSize>1)
|
|
|
|
{
|
2011-02-24 05:22:10 +08:00
|
|
|
for(int offset=0;offset<4;++offset)
|
|
|
|
{
|
|
|
|
for(int i=0;i<PacketSize/2;++i)
|
|
|
|
ref[2*i+0] = ref[2*i+1] = data1[offset+i];
|
|
|
|
internal::pstore(data2,internal::ploaddup<Packet>(data1+offset));
|
|
|
|
VERIFY(areApprox(ref, data2, PacketSize) && "ploaddup");
|
|
|
|
}
|
2011-02-23 23:20:55 +08:00
|
|
|
}
|
2008-08-21 04:08:38 +08:00
|
|
|
|
|
|
|
ref[0] = 0;
|
|
|
|
for (int i=0; i<PacketSize; ++i)
|
|
|
|
ref[0] += data1[i];
|
2010-10-25 22:15:22 +08:00
|
|
|
VERIFY(isApproxAbs(ref[0], internal::predux(internal::pload<Packet>(data1)), refvalue) && "internal::predux");
|
2009-03-10 02:40:09 +08:00
|
|
|
|
2009-02-11 02:06:05 +08:00
|
|
|
ref[0] = 1;
|
|
|
|
for (int i=0; i<PacketSize; ++i)
|
|
|
|
ref[0] *= data1[i];
|
2010-10-25 22:15:22 +08:00
|
|
|
VERIFY(internal::isApprox(ref[0], internal::predux_mul(internal::pload<Packet>(data1))) && "internal::predux_mul");
|
2009-03-10 02:40:09 +08:00
|
|
|
|
2008-08-21 04:08:38 +08:00
|
|
|
for (int j=0; j<PacketSize; ++j)
|
|
|
|
{
|
|
|
|
ref[j] = 0;
|
|
|
|
for (int i=0; i<PacketSize; ++i)
|
|
|
|
ref[j] += data1[i+j*PacketSize];
|
2010-10-25 22:15:22 +08:00
|
|
|
packets[j] = internal::pload<Packet>(data1+j*PacketSize);
|
2008-08-21 04:08:38 +08:00
|
|
|
}
|
2010-10-25 22:15:22 +08:00
|
|
|
internal::pstore(data2, internal::preduxp(packets));
|
|
|
|
VERIFY(areApproxAbs(ref, data2, PacketSize, refvalue) && "internal::preduxp");
|
2009-02-06 20:40:38 +08:00
|
|
|
|
|
|
|
for (int i=0; i<PacketSize; ++i)
|
|
|
|
ref[i] = data1[PacketSize-i-1];
|
2010-10-25 22:15:22 +08:00
|
|
|
internal::pstore(data2, internal::preverse(internal::pload<Packet>(data1)));
|
|
|
|
VERIFY(areApprox(ref, data2, PacketSize) && "internal::preverse");
|
2008-08-21 04:08:38 +08:00
|
|
|
}
|
|
|
|
|
2009-03-25 20:26:13 +08:00
|
|
|
template<typename Scalar> void packetmath_real()
|
|
|
|
{
|
2012-11-06 22:25:50 +08:00
|
|
|
using std::abs;
|
2010-10-25 22:15:22 +08:00
|
|
|
typedef typename internal::packet_traits<Scalar>::type Packet;
|
|
|
|
const int PacketSize = internal::packet_traits<Scalar>::size;
|
2009-03-25 20:26:13 +08:00
|
|
|
|
|
|
|
const int size = PacketSize*4;
|
2010-10-25 22:15:22 +08:00
|
|
|
EIGEN_ALIGN16 Scalar data1[internal::packet_traits<Scalar>::size*4];
|
|
|
|
EIGEN_ALIGN16 Scalar data2[internal::packet_traits<Scalar>::size*4];
|
|
|
|
EIGEN_ALIGN16 Scalar ref[internal::packet_traits<Scalar>::size*4];
|
2010-07-05 16:54:24 +08:00
|
|
|
|
2009-03-25 20:26:13 +08:00
|
|
|
for (int i=0; i<size; ++i)
|
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
data1[i] = internal::random<Scalar>(-1e3,1e3);
|
|
|
|
data2[i] = internal::random<Scalar>(-1e3,1e3);
|
2009-03-25 20:26:13 +08:00
|
|
|
}
|
2012-11-06 22:25:50 +08:00
|
|
|
CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasSin, std::sin, internal::psin);
|
|
|
|
CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasCos, std::cos, internal::pcos);
|
|
|
|
CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasTan, std::tan, internal::ptan);
|
2011-02-18 00:37:11 +08:00
|
|
|
|
|
|
|
for (int i=0; i<size; ++i)
|
|
|
|
{
|
|
|
|
data1[i] = internal::random<Scalar>(-1,1);
|
|
|
|
data2[i] = internal::random<Scalar>(-1,1);
|
|
|
|
}
|
2012-11-06 22:25:50 +08:00
|
|
|
CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasASin, std::asin, internal::pasin);
|
|
|
|
CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasACos, std::acos, internal::pacos);
|
2010-07-05 16:54:24 +08:00
|
|
|
|
2009-03-25 20:26:13 +08:00
|
|
|
for (int i=0; i<size; ++i)
|
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
data1[i] = internal::random<Scalar>(-87,88);
|
|
|
|
data2[i] = internal::random<Scalar>(-87,88);
|
2009-03-25 20:26:13 +08:00
|
|
|
}
|
2012-11-06 22:25:50 +08:00
|
|
|
CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasExp, std::exp, internal::pexp);
|
2010-07-05 16:54:24 +08:00
|
|
|
|
2009-03-25 20:26:13 +08:00
|
|
|
for (int i=0; i<size; ++i)
|
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
data1[i] = internal::random<Scalar>(0,1e6);
|
|
|
|
data2[i] = internal::random<Scalar>(0,1e6);
|
2009-03-25 20:26:13 +08:00
|
|
|
}
|
2013-02-15 06:34:05 +08:00
|
|
|
if(internal::random<float>(0,1)<0.1)
|
|
|
|
data1[internal::random<int>(0, PacketSize)] = 0;
|
2012-11-06 22:25:50 +08:00
|
|
|
CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasLog, std::log, internal::plog);
|
|
|
|
CHECK_CWISE1_IF(internal::packet_traits<Scalar>::HasSqrt, std::sqrt, internal::psqrt);
|
2013-03-21 01:28:40 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
template<typename Scalar> void packetmath_notcomplex()
|
|
|
|
{
|
|
|
|
using std::abs;
|
|
|
|
typedef typename internal::packet_traits<Scalar>::type Packet;
|
|
|
|
const int PacketSize = internal::packet_traits<Scalar>::size;
|
|
|
|
|
|
|
|
EIGEN_ALIGN16 Scalar data1[internal::packet_traits<Scalar>::size*4];
|
|
|
|
EIGEN_ALIGN16 Scalar data2[internal::packet_traits<Scalar>::size*4];
|
|
|
|
EIGEN_ALIGN16 Scalar ref[internal::packet_traits<Scalar>::size*4];
|
2013-04-10 15:46:16 +08:00
|
|
|
|
|
|
|
Array<Scalar,Dynamic,1>::Map(data1, internal::packet_traits<Scalar>::size*4).setRandom();
|
2010-07-05 22:18:09 +08:00
|
|
|
|
|
|
|
ref[0] = data1[0];
|
|
|
|
for (int i=0; i<PacketSize; ++i)
|
2011-08-19 20:18:05 +08:00
|
|
|
ref[0] = (std::min)(ref[0],data1[i]);
|
2010-10-25 22:15:22 +08:00
|
|
|
VERIFY(internal::isApprox(ref[0], internal::predux_min(internal::pload<Packet>(data1))) && "internal::predux_min");
|
2010-07-05 22:18:09 +08:00
|
|
|
|
2011-08-19 20:18:05 +08:00
|
|
|
CHECK_CWISE2((std::min), internal::pmin);
|
|
|
|
CHECK_CWISE2((std::max), internal::pmax);
|
2012-11-06 22:25:50 +08:00
|
|
|
CHECK_CWISE1(abs, internal::pabs);
|
2010-07-05 22:18:09 +08:00
|
|
|
|
|
|
|
ref[0] = data1[0];
|
|
|
|
for (int i=0; i<PacketSize; ++i)
|
2011-08-19 20:18:05 +08:00
|
|
|
ref[0] = (std::max)(ref[0],data1[i]);
|
2010-10-25 22:15:22 +08:00
|
|
|
VERIFY(internal::isApprox(ref[0], internal::predux_max(internal::pload<Packet>(data1))) && "internal::predux_max");
|
2011-05-19 03:11:03 +08:00
|
|
|
|
|
|
|
for (int i=0; i<PacketSize; ++i)
|
|
|
|
ref[i] = data1[0]+Scalar(i);
|
|
|
|
internal::pstore(data2, internal::plset(data1[0]));
|
|
|
|
VERIFY(areApprox(ref, data2, PacketSize) && "internal::plset");
|
2009-03-25 20:26:13 +08:00
|
|
|
}
|
|
|
|
|
2011-02-24 02:24:26 +08:00
|
|
|
template<typename Scalar,bool ConjLhs,bool ConjRhs> void test_conj_helper(Scalar* data1, Scalar* data2, Scalar* ref, Scalar* pval)
|
|
|
|
{
|
|
|
|
typedef typename internal::packet_traits<Scalar>::type Packet;
|
|
|
|
const int PacketSize = internal::packet_traits<Scalar>::size;
|
|
|
|
|
|
|
|
internal::conj_if<ConjLhs> cj0;
|
|
|
|
internal::conj_if<ConjRhs> cj1;
|
|
|
|
internal::conj_helper<Scalar,Scalar,ConjLhs,ConjRhs> cj;
|
|
|
|
internal::conj_helper<Packet,Packet,ConjLhs,ConjRhs> pcj;
|
|
|
|
|
|
|
|
for(int i=0;i<PacketSize;++i)
|
|
|
|
{
|
|
|
|
ref[i] = cj0(data1[i]) * cj1(data2[i]);
|
|
|
|
VERIFY(internal::isApprox(ref[i], cj.pmul(data1[i],data2[i])) && "conj_helper pmul");
|
|
|
|
}
|
|
|
|
internal::pstore(pval,pcj.pmul(internal::pload<Packet>(data1),internal::pload<Packet>(data2)));
|
|
|
|
VERIFY(areApprox(ref, pval, PacketSize) && "conj_helper pmul");
|
|
|
|
|
|
|
|
for(int i=0;i<PacketSize;++i)
|
|
|
|
{
|
|
|
|
Scalar tmp = ref[i];
|
|
|
|
ref[i] += cj0(data1[i]) * cj1(data2[i]);
|
|
|
|
VERIFY(internal::isApprox(ref[i], cj.pmadd(data1[i],data2[i],tmp)) && "conj_helper pmadd");
|
|
|
|
}
|
|
|
|
internal::pstore(pval,pcj.pmadd(internal::pload<Packet>(data1),internal::pload<Packet>(data2),internal::pload<Packet>(pval)));
|
|
|
|
VERIFY(areApprox(ref, pval, PacketSize) && "conj_helper pmadd");
|
|
|
|
}
|
|
|
|
|
2010-07-07 02:54:14 +08:00
|
|
|
template<typename Scalar> void packetmath_complex()
|
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
typedef typename internal::packet_traits<Scalar>::type Packet;
|
|
|
|
const int PacketSize = internal::packet_traits<Scalar>::size;
|
2010-07-07 02:54:14 +08:00
|
|
|
|
|
|
|
const int size = PacketSize*4;
|
|
|
|
EIGEN_ALIGN16 Scalar data1[PacketSize*4];
|
|
|
|
EIGEN_ALIGN16 Scalar data2[PacketSize*4];
|
|
|
|
EIGEN_ALIGN16 Scalar ref[PacketSize*4];
|
|
|
|
EIGEN_ALIGN16 Scalar pval[PacketSize*4];
|
|
|
|
|
|
|
|
for (int i=0; i<size; ++i)
|
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
data1[i] = internal::random<Scalar>() * Scalar(1e2);
|
|
|
|
data2[i] = internal::random<Scalar>() * Scalar(1e2);
|
2010-07-07 02:54:14 +08:00
|
|
|
}
|
2011-02-24 02:24:26 +08:00
|
|
|
|
|
|
|
test_conj_helper<Scalar,false,false> (data1,data2,ref,pval);
|
|
|
|
test_conj_helper<Scalar,false,true> (data1,data2,ref,pval);
|
|
|
|
test_conj_helper<Scalar,true,false> (data1,data2,ref,pval);
|
|
|
|
test_conj_helper<Scalar,true,true> (data1,data2,ref,pval);
|
2010-07-07 02:54:14 +08:00
|
|
|
|
2011-02-23 21:20:33 +08:00
|
|
|
{
|
|
|
|
for(int i=0;i<PacketSize;++i)
|
|
|
|
ref[i] = Scalar(std::imag(data1[i]),std::real(data1[i]));
|
|
|
|
internal::pstore(pval,internal::pcplxflip(internal::pload<Packet>(data1)));
|
|
|
|
VERIFY(areApprox(ref, pval, PacketSize) && "pcplxflip");
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2010-07-07 02:54:14 +08:00
|
|
|
}
|
|
|
|
|
2008-08-21 04:08:38 +08:00
|
|
|
void test_packetmath()
|
|
|
|
{
|
|
|
|
for(int i = 0; i < g_repeat; i++) {
|
2010-07-07 02:54:14 +08:00
|
|
|
CALL_SUBTEST_1( packetmath<float>() );
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_2( packetmath<double>() );
|
|
|
|
CALL_SUBTEST_3( packetmath<int>() );
|
|
|
|
CALL_SUBTEST_1( packetmath<std::complex<float> >() );
|
2010-07-07 02:54:14 +08:00
|
|
|
CALL_SUBTEST_2( packetmath<std::complex<double> >() );
|
2010-07-05 16:54:24 +08:00
|
|
|
|
2013-03-21 01:28:40 +08:00
|
|
|
CALL_SUBTEST_1( packetmath_notcomplex<float>() );
|
|
|
|
CALL_SUBTEST_2( packetmath_notcomplex<double>() );
|
|
|
|
CALL_SUBTEST_3( packetmath_notcomplex<int>() );
|
|
|
|
|
2010-07-07 02:54:14 +08:00
|
|
|
CALL_SUBTEST_1( packetmath_real<float>() );
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_2( packetmath_real<double>() );
|
2010-07-07 02:54:14 +08:00
|
|
|
|
|
|
|
CALL_SUBTEST_1( packetmath_complex<std::complex<float> >() );
|
|
|
|
CALL_SUBTEST_2( packetmath_complex<std::complex<double> >() );
|
2008-08-21 04:08:38 +08:00
|
|
|
}
|
|
|
|
}
|