aboutsummaryrefslogtreecommitdiffhomepage
path: root/test
diff options
context:
space:
mode:
authorGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2009-10-15 16:09:43 -0400
committerGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2009-10-15 16:09:43 -0400
commit3c4a025a5461262e175542fa8857e2bb6e0c63b2 (patch)
tree6cd201cb666b53c9c7817d0a00c21569bc4674ab /test
parent41e942d3fb3a3bf26a7fad169adaf3968daa7d46 (diff)
parentd177c1f3aca150d36d17c8116504e3ea1e7b30cf (diff)
merge
Diffstat (limited to 'test')
-rw-r--r--test/CMakeLists.txt4
-rw-r--r--test/array_replicate.cpp10
-rw-r--r--test/main.h12
-rw-r--r--test/packetmath.cpp14
-rw-r--r--test/product_extra.cpp1
-rw-r--r--test/product_small.cpp6
-rw-r--r--test/qr.cpp44
-rw-r--r--test/qr_colpivoting.cpp42
-rw-r--r--test/qr_fullpivoting.cpp6
-rw-r--r--test/stable_norm.cpp14
-rw-r--r--test/testsuite.cmake6
-rw-r--r--test/unalignedassert.cpp96
-rw-r--r--test/visitor.cpp4
13 files changed, 173 insertions, 86 deletions
diff --git a/test/CMakeLists.txt b/test/CMakeLists.txt
index f3c15612f..b46971493 100644
--- a/test/CMakeLists.txt
+++ b/test/CMakeLists.txt
@@ -1,6 +1,6 @@
-
+project(EigenTesting)
+add_custom_target(btest)
include(EigenTesting)
-enable_testing()
ei_init_testing()
find_package(GSL)
diff --git a/test/array_replicate.cpp b/test/array_replicate.cpp
index d1608915f..cd0f65f26 100644
--- a/test/array_replicate.cpp
+++ b/test/array_replicate.cpp
@@ -42,9 +42,9 @@ template<typename MatrixType> void replicate(const MatrixType& m)
MatrixType m1 = MatrixType::Random(rows, cols),
m2 = MatrixType::Random(rows, cols);
-
+
VectorType v1 = VectorType::Random(rows);
-
+
MatrixX x1, x2;
VectorX vx1;
@@ -56,17 +56,17 @@ template<typename MatrixType> void replicate(const MatrixType& m)
for(int i=0; i<f1; i++)
x1.block(i*rows,j*cols,rows,cols) = m1;
VERIFY_IS_APPROX(x1, m1.replicate(f1,f2));
-
+
x2.resize(2*rows,3*cols);
x2 << m2, m2, m2,
m2, m2, m2;
VERIFY_IS_APPROX(x2, (m2.template replicate<2,3>()));
-
+
x2.resize(rows,f1);
for (int j=0; j<f1; ++j)
x2.col(j) = v1;
VERIFY_IS_APPROX(x2, v1.rowwise().replicate(f1));
-
+
vx1.resize(rows*f2);
for (int j=0; j<f2; ++j)
vx1.segment(j*rows,rows) = v1;
diff --git a/test/main.h b/test/main.h
index 8c93e856c..51b719814 100644
--- a/test/main.h
+++ b/test/main.h
@@ -248,18 +248,22 @@ template<typename MatrixType>
void createRandomMatrixOfRank(int desired_rank, int rows, int cols, MatrixType& m)
{
typedef typename ei_traits<MatrixType>::Scalar Scalar;
+ enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime };
+
typedef Matrix<Scalar, Dynamic, 1> VectorType;
+ typedef Matrix<Scalar, Rows, Rows> MatrixAType;
+ typedef Matrix<Scalar, Cols, Cols> MatrixBType;
- MatrixType a = MatrixType::Random(rows,rows);
+ MatrixAType a = MatrixAType::Random(rows,rows);
MatrixType d = MatrixType::Identity(rows,cols);
- MatrixType b = MatrixType::Random(cols,cols);
+ MatrixBType b = MatrixBType::Random(cols,cols);
// set the diagonal such that only desired_rank non-zero entries reamain
const int diag_size = std::min(d.rows(),d.cols());
d.diagonal().segment(desired_rank, diag_size-desired_rank) = VectorType::Zero(diag_size-desired_rank);
- HouseholderQR<MatrixType> qra(a);
- HouseholderQR<MatrixType> qrb(b);
+ HouseholderQR<MatrixAType> qra(a);
+ HouseholderQR<MatrixBType> qrb(b);
m = qra.matrixQ() * d * qrb.matrixQ();
}
diff --git a/test/packetmath.cpp b/test/packetmath.cpp
index d86d40d68..1745ae5c6 100644
--- a/test/packetmath.cpp
+++ b/test/packetmath.cpp
@@ -99,10 +99,10 @@ template<typename Scalar> void packetmath()
const int PacketSize = ei_packet_traits<Scalar>::size;
const int size = PacketSize*4;
- EIGEN_ALIGN_128 Scalar data1[ei_packet_traits<Scalar>::size*4];
- EIGEN_ALIGN_128 Scalar data2[ei_packet_traits<Scalar>::size*4];
- EIGEN_ALIGN_128 Packet packets[PacketSize*2];
- EIGEN_ALIGN_128 Scalar ref[ei_packet_traits<Scalar>::size*4];
+ EIGEN_ALIGN16 Scalar data1[ei_packet_traits<Scalar>::size*4];
+ EIGEN_ALIGN16 Scalar data2[ei_packet_traits<Scalar>::size*4];
+ EIGEN_ALIGN16 Packet packets[PacketSize*2];
+ EIGEN_ALIGN16 Scalar ref[ei_packet_traits<Scalar>::size*4];
for (int i=0; i<size; ++i)
{
data1[i] = ei_random<Scalar>();
@@ -202,9 +202,9 @@ template<typename Scalar> void packetmath_real()
const int PacketSize = ei_packet_traits<Scalar>::size;
const int size = PacketSize*4;
- EIGEN_ALIGN_128 Scalar data1[ei_packet_traits<Scalar>::size*4];
- EIGEN_ALIGN_128 Scalar data2[ei_packet_traits<Scalar>::size*4];
- EIGEN_ALIGN_128 Scalar ref[ei_packet_traits<Scalar>::size*4];
+ EIGEN_ALIGN16 Scalar data1[ei_packet_traits<Scalar>::size*4];
+ EIGEN_ALIGN16 Scalar data2[ei_packet_traits<Scalar>::size*4];
+ EIGEN_ALIGN16 Scalar ref[ei_packet_traits<Scalar>::size*4];
for (int i=0; i<size; ++i)
{
diff --git a/test/product_extra.cpp b/test/product_extra.cpp
index 3ad99fc7a..8e55c6010 100644
--- a/test/product_extra.cpp
+++ b/test/product_extra.cpp
@@ -114,7 +114,6 @@ template<typename MatrixType> void product_extra(const MatrixType& m)
VERIFY_IS_APPROX(m1.col(j2).adjoint() * m1.block(0,j,m1.rows(),c), m1.col(j2).adjoint().eval() * m1.block(0,j,m1.rows(),c).eval());
VERIFY_IS_APPROX(m1.block(i,0,r,m1.cols()) * m1.row(i2).adjoint(), m1.block(i,0,r,m1.cols()).eval() * m1.row(i2).adjoint().eval());
-
}
void test_product_extra()
diff --git a/test/product_small.cpp b/test/product_small.cpp
index 3aed5cf1b..182af71db 100644
--- a/test/product_small.cpp
+++ b/test/product_small.cpp
@@ -34,4 +34,10 @@ void test_product_small()
CALL_SUBTEST( product(Matrix4d()) );
CALL_SUBTEST( product(Matrix4f()) );
}
+
+ {
+ // test compilation of (outer_product) * vector
+ Vector3f v = Vector3f::Random();
+ VERIFY_IS_APPROX( (v * v.transpose()) * v, (v * v.transpose()).eval() * v);
+ }
}
diff --git a/test/qr.cpp b/test/qr.cpp
index f185ac86e..864828750 100644
--- a/test/qr.cpp
+++ b/test/qr.cpp
@@ -31,30 +31,40 @@ template<typename MatrixType> void qr(const MatrixType& m)
int cols = m.cols();
typedef typename MatrixType::Scalar Scalar;
- typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> SquareMatrixType;
+ typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
MatrixType a = MatrixType::Random(rows,cols);
HouseholderQR<MatrixType> qrOfA(a);
MatrixType r = qrOfA.matrixQR();
+
+ MatrixQType q = qrOfA.matrixQ();
+ VERIFY_IS_UNITARY(q);
+
// FIXME need better way to construct trapezoid
for(int i = 0; i < rows; i++) for(int j = 0; j < cols; j++) if(i>j) r(i,j) = Scalar(0);
VERIFY_IS_APPROX(a, qrOfA.matrixQ() * r);
+}
- SquareMatrixType b = a.adjoint() * a;
+template<typename MatrixType, int Cols2> void qr_fixedsize()
+{
+ enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime };
+ typedef typename MatrixType::Scalar Scalar;
+ Matrix<Scalar,Rows,Cols> m1 = Matrix<Scalar,Rows,Cols>::Random();
+ HouseholderQR<Matrix<Scalar,Rows,Cols> > qr(m1);
- // check tridiagonalization
- Tridiagonalization<SquareMatrixType> tridiag(b);
- VERIFY_IS_APPROX(b, tridiag.matrixQ() * tridiag.matrixT() * tridiag.matrixQ().adjoint());
+ Matrix<Scalar,Rows,Cols> r = qr.matrixQR();
+ // FIXME need better way to construct trapezoid
+ for(int i = 0; i < Rows; i++) for(int j = 0; j < Cols; j++) if(i>j) r(i,j) = Scalar(0);
- // check hessenberg decomposition
- HessenbergDecomposition<SquareMatrixType> hess(b);
- VERIFY_IS_APPROX(b, hess.matrixQ() * hess.matrixH() * hess.matrixQ().adjoint());
- VERIFY_IS_APPROX(tridiag.matrixT(), hess.matrixH());
- b = SquareMatrixType::Random(cols,cols);
- hess.compute(b);
- VERIFY_IS_APPROX(b, hess.matrixQ() * hess.matrixH() * hess.matrixQ().adjoint());
+ VERIFY_IS_APPROX(m1, qr.matrixQ() * r);
+
+ Matrix<Scalar,Cols,Cols2> m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2);
+ Matrix<Scalar,Rows,Cols2> m3 = m1*m2;
+ m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2);
+ qr.solve(m3, &m2);
+ VERIFY_IS_APPROX(m3, m1*m2);
}
template<typename MatrixType> void qr_invertible()
@@ -105,11 +115,11 @@ template<typename MatrixType> void qr_verify_assert()
void test_qr()
{
for(int i = 0; i < 1; i++) {
- // FIXME : very weird bug here
-// CALL_SUBTEST( qr(Matrix2f()) );
- CALL_SUBTEST( qr(Matrix4d()) );
- CALL_SUBTEST( qr(MatrixXf(47,40)) );
- CALL_SUBTEST( qr(MatrixXcd(17,7)) );
+ CALL_SUBTEST( qr(MatrixXf(47,40)) );
+ CALL_SUBTEST( qr(MatrixXcd(17,7)) );
+ CALL_SUBTEST(( qr_fixedsize<Matrix<float,3,4>, 2 >() ));
+ CALL_SUBTEST(( qr_fixedsize<Matrix<double,6,2>, 4 >() ));
+ CALL_SUBTEST(( qr_fixedsize<Matrix<double,2,5>, 7 >() ));
}
for(int i = 0; i < g_repeat; i++) {
diff --git a/test/qr_colpivoting.cpp b/test/qr_colpivoting.cpp
index 588a41e56..5c5c5d259 100644
--- a/test/qr_colpivoting.cpp
+++ b/test/qr_colpivoting.cpp
@@ -32,7 +32,7 @@ template<typename MatrixType> void qr()
int rank = ei_random<int>(1, std::min(rows, cols)-1);
typedef typename MatrixType::Scalar Scalar;
- typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> SquareMatrixType;
+ typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
MatrixType m1;
createRandomMatrixOfRank(rank,rows,cols,m1);
@@ -44,6 +44,10 @@ template<typename MatrixType> void qr()
VERIFY(!qr.isSurjective());
MatrixType r = qr.matrixQR();
+
+ MatrixQType q = qr.matrixQ();
+ VERIFY_IS_UNITARY(q);
+
// FIXME need better way to construct trapezoid
for(int i = 0; i < rows; i++) for(int j = 0; j < cols; j++) if(i>j) r(i,j) = Scalar(0);
@@ -63,6 +67,40 @@ template<typename MatrixType> void qr()
VERIFY(!qr.solve(m3, &m2));
}
+template<typename MatrixType, int Cols2> void qr_fixedsize()
+{
+ enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime };
+ typedef typename MatrixType::Scalar Scalar;
+ int rank = ei_random<int>(1, std::min(int(Rows), int(Cols))-1);
+ Matrix<Scalar,Rows,Cols> m1;
+ createRandomMatrixOfRank(rank,Rows,Cols,m1);
+ ColPivotingHouseholderQR<Matrix<Scalar,Rows,Cols> > qr(m1);
+ VERIFY_IS_APPROX(rank, qr.rank());
+ VERIFY(Cols - qr.rank() == qr.dimensionOfKernel());
+ VERIFY(!qr.isInjective());
+ VERIFY(!qr.isInvertible());
+ VERIFY(!qr.isSurjective());
+
+ Matrix<Scalar,Rows,Cols> r = qr.matrixQR();
+ // FIXME need better way to construct trapezoid
+ for(int i = 0; i < Rows; i++) for(int j = 0; j < Cols; j++) if(i>j) r(i,j) = Scalar(0);
+
+ Matrix<Scalar,Rows,Cols> b = qr.matrixQ() * r;
+
+ Matrix<Scalar,Rows,Cols> c = MatrixType::Zero(Rows,Cols);
+
+ for(int i = 0; i < Cols; ++i) c.col(qr.colsPermutation().coeff(i)) = b.col(i);
+ VERIFY_IS_APPROX(m1, c);
+
+ Matrix<Scalar,Cols,Cols2> m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2);
+ Matrix<Scalar,Rows,Cols2> m3 = m1*m2;
+ m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2);
+ VERIFY(qr.solve(m3, &m2));
+ VERIFY_IS_APPROX(m3, m1*m2);
+ m3 = Matrix<Scalar,Rows,Cols2>::Random(Rows,Cols2);
+ VERIFY(!qr.solve(m3, &m2));
+}
+
template<typename MatrixType> void qr_invertible()
{
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
@@ -120,6 +158,8 @@ void test_qr_colpivoting()
CALL_SUBTEST( qr<MatrixXf>() );
CALL_SUBTEST( qr<MatrixXd>() );
CALL_SUBTEST( qr<MatrixXcd>() );
+ CALL_SUBTEST(( qr_fixedsize<Matrix<float,3,5>, 4 >() ));
+ CALL_SUBTEST(( qr_fixedsize<Matrix<double,6,2>, 3 >() ));
}
for(int i = 0; i < g_repeat; i++) {
diff --git a/test/qr_fullpivoting.cpp b/test/qr_fullpivoting.cpp
index 3a37bcb46..891c2a527 100644
--- a/test/qr_fullpivoting.cpp
+++ b/test/qr_fullpivoting.cpp
@@ -32,7 +32,7 @@ template<typename MatrixType> void qr()
int rank = ei_random<int>(1, std::min(rows, cols)-1);
typedef typename MatrixType::Scalar Scalar;
- typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> SquareMatrixType;
+ typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
MatrixType m1;
createRandomMatrixOfRank(rank,rows,cols,m1);
@@ -44,6 +44,10 @@ template<typename MatrixType> void qr()
VERIFY(!qr.isSurjective());
MatrixType r = qr.matrixQR();
+
+ MatrixQType q = qr.matrixQ();
+ VERIFY_IS_UNITARY(q);
+
// FIXME need better way to construct trapezoid
for(int i = 0; i < rows; i++) for(int j = 0; j < cols; j++) if(i>j) r(i,j) = Scalar(0);
diff --git a/test/stable_norm.cpp b/test/stable_norm.cpp
index 726512ec0..ed72bb7a7 100644
--- a/test/stable_norm.cpp
+++ b/test/stable_norm.cpp
@@ -33,6 +33,20 @@ template<typename MatrixType> void stable_norm(const MatrixType& m)
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
+ // Check the basic machine-dependent constants.
+ {
+ int ibeta, it, iemin, iemax;
+
+ ibeta = std::numeric_limits<RealScalar>::radix; // base for floating-point numbers
+ it = std::numeric_limits<RealScalar>::digits; // number of base-beta digits in mantissa
+ iemin = std::numeric_limits<RealScalar>::min_exponent; // minimum exponent
+ iemax = std::numeric_limits<RealScalar>::max_exponent; // maximum exponent
+
+ VERIFY( (!(iemin > 1 - 2*it || 1+it>iemax || (it==2 && ibeta<5) || (it<=4 && ibeta <= 3 ) || it<2))
+ && "the stable norm algorithm cannot be guaranteed on this computer");
+ }
+
+
int rows = m.rows();
int cols = m.cols();
diff --git a/test/testsuite.cmake b/test/testsuite.cmake
index 37ee87565..5d0cb6585 100644
--- a/test/testsuite.cmake
+++ b/test/testsuite.cmake
@@ -148,7 +148,7 @@ endif(NOT EIGEN_NO_UPDATE)
# which ctest command to use for running the dashboard
SET (CTEST_COMMAND "${EIGEN_CMAKE_DIR}ctest -D ${EIGEN_MODE}")
# what cmake command to use for configuring this dashboard
-SET (CTEST_CMAKE_COMMAND "${EIGEN_CMAKE_DIR}cmake -DEIGEN_BUILD_TESTS=on ")
+SET (CTEST_CMAKE_COMMAND "${EIGEN_CMAKE_DIR}cmake -DEIGEN_CMAKE_RUN_FROM_CTEST=ON")
####################################################################
# The values in this section are optional you can either
@@ -175,9 +175,9 @@ if(WIN32 AND NOT UNIX)
else(EIGEN_GENERATOR_TYPE)
set(CTEST_CMAKE_COMMAND "${CTEST_CMAKE_COMMAND} -G \"NMake Makefiles\" -DCMAKE_MAKE_PROGRAM=nmake")
SET (CTEST_INITIAL_CACHE "
- MAKECOMMAND:STRING=nmake -i
+ MAKECOMMAND:STRING=nmake /i
CMAKE_MAKE_PROGRAM:FILEPATH=nmake
- CMAKE_GENERATOR:INTERNAL=NMake Makefiles
+ CMAKE_GENERATOR:INTERNAL=NMake Makefiles
CMAKE_BUILD_TYPE:STRING=Release
BUILDNAME:STRING=${EIGEN_BUILD_STRING}
SITE:STRING=${EIGEN_SITE}
diff --git a/test/unalignedassert.cpp b/test/unalignedassert.cpp
index ade1ab26e..2b819417e 100644
--- a/test/unalignedassert.cpp
+++ b/test/unalignedassert.cpp
@@ -24,52 +24,38 @@
#include "main.h"
-struct Good1
+struct TestNew1
{
MatrixXd m; // good: m will allocate its own array, taking care of alignment.
- Good1() : m(20,20) {}
+ TestNew1() : m(20,20) {}
};
-struct Good2
+struct TestNew2
{
- Matrix3d m; // good: m's size isn't a multiple of 16 bytes, so m doesn't have to be aligned
+ Matrix3d m; // good: m's size isn't a multiple of 16 bytes, so m doesn't have to be 16-byte aligned,
+ // 8-byte alignment is good enough here, which we'll get automatically
};
-struct Good3
+struct TestNew3
{
- Vector2f m; // good: same reason
+ Vector2f m; // good: m's size isn't a multiple of 16 bytes, so m doesn't have to be 16-byte aligned
};
-struct Bad4
-{
- Vector2d m; // bad: sizeof(m)%16==0 so alignment is required
-};
-
-struct Bad5
-{
- Matrix<float, 2, 6> m; // bad: same reason
-};
-
-struct Bad6
-{
- Matrix<double, 3, 4> m; // bad: same reason
-};
-
-struct Good7
+struct TestNew4
{
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
Vector2d m;
float f; // make the struct have sizeof%16!=0 to make it a little more tricky when we allow an array of 2 such objects
};
-struct Good8
+struct TestNew5
{
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
- float f; // try the f at first -- the EIGEN_ALIGN_128 attribute of m should make that still work
+ float f; // try the f at first -- the EIGEN_ALIGN16 attribute of m should make that still work
Matrix4f m;
};
-struct Good9
+struct TestNew6
{
Matrix<float,2,2,DontAlign> m; // good: no alignment requested
float f;
@@ -94,34 +80,58 @@ void check_unalignedassert_good()
#if EIGEN_ALIGN
template<typename T>
-void check_unalignedassert_bad()
+void construct_at_boundary(int boundary)
{
- float buf[sizeof(T)+16];
- float *unaligned = buf;
- while((reinterpret_cast<size_t>(unaligned)&0xf)==0) ++unaligned; // make sure unaligned is really unaligned
- T *x = ::new(static_cast<void*>(unaligned)) T;
+ char buf[sizeof(T)+256];
+ size_t _buf = reinterpret_cast<size_t>(buf);
+ _buf += (16 - (_buf % 16)); // make 16-byte aligned
+ _buf += boundary; // make exact boundary-aligned
+ T *x = ::new(reinterpret_cast<void*>(_buf)) T;
x->~T();
}
#endif
void unalignedassert()
{
- check_unalignedassert_good<Good1>();
- check_unalignedassert_good<Good2>();
- check_unalignedassert_good<Good3>();
-#if EIGEN_ALIGN
- VERIFY_RAISES_ASSERT(check_unalignedassert_bad<Bad4>());
- VERIFY_RAISES_ASSERT(check_unalignedassert_bad<Bad5>());
- VERIFY_RAISES_ASSERT(check_unalignedassert_bad<Bad6>());
-#endif
-
- check_unalignedassert_good<Good7>();
- check_unalignedassert_good<Good8>();
- check_unalignedassert_good<Good9>();
+ #if EIGEN_ALIGN
+ construct_at_boundary<Vector2f>(4);
+ construct_at_boundary<Vector3f>(4);
+ construct_at_boundary<Vector4f>(16);
+ construct_at_boundary<Matrix2f>(16);
+ construct_at_boundary<Matrix3f>(4);
+ construct_at_boundary<Matrix4f>(16);
+
+ construct_at_boundary<Vector2d>(16);
+ construct_at_boundary<Vector3d>(4);
+ construct_at_boundary<Vector4d>(16);
+ construct_at_boundary<Matrix2d>(16);
+ construct_at_boundary<Matrix3d>(4);
+ construct_at_boundary<Matrix4d>(16);
+
+ construct_at_boundary<Vector2cf>(16);
+ construct_at_boundary<Vector3cf>(4);
+ construct_at_boundary<Vector2cd>(16);
+ construct_at_boundary<Vector3cd>(16);
+ #endif
+
+ check_unalignedassert_good<TestNew1>();
+ check_unalignedassert_good<TestNew2>();
+ check_unalignedassert_good<TestNew3>();
+
+ check_unalignedassert_good<TestNew4>();
+ check_unalignedassert_good<TestNew5>();
+ check_unalignedassert_good<TestNew6>();
check_unalignedassert_good<Depends<true> >();
#if EIGEN_ALIGN
- VERIFY_RAISES_ASSERT(check_unalignedassert_bad<Depends<false> >());
+ VERIFY_RAISES_ASSERT(construct_at_boundary<Vector4f>(8));
+ VERIFY_RAISES_ASSERT(construct_at_boundary<Matrix4f>(8));
+ VERIFY_RAISES_ASSERT(construct_at_boundary<Vector2d>(8));
+ VERIFY_RAISES_ASSERT(construct_at_boundary<Vector4d>(8));
+ VERIFY_RAISES_ASSERT(construct_at_boundary<Matrix2d>(8));
+ VERIFY_RAISES_ASSERT(construct_at_boundary<Matrix4d>(8));
+ VERIFY_RAISES_ASSERT(construct_at_boundary<Vector2cf>(8));
+ VERIFY_RAISES_ASSERT(construct_at_boundary<Vector2cd>(8));
#endif
}
diff --git a/test/visitor.cpp b/test/visitor.cpp
index b78782b78..6ec442bc8 100644
--- a/test/visitor.cpp
+++ b/test/visitor.cpp
@@ -40,7 +40,7 @@ template<typename MatrixType> void matrixVisitor(const MatrixType& p)
m(i) = ei_random<Scalar>();
Scalar minc = Scalar(1000), maxc = Scalar(-1000);
- int minrow,mincol,maxrow,maxcol;
+ int minrow=0,mincol=0,maxrow=0,maxcol=0;
for(int j = 0; j < cols; j++)
for(int i = 0; i < rows; i++)
{
@@ -86,7 +86,7 @@ template<typename VectorType> void vectorVisitor(const VectorType& w)
v(i) = ei_random<Scalar>();
Scalar minc = Scalar(1000), maxc = Scalar(-1000);
- int minidx,maxidx;
+ int minidx=0,maxidx=0;
for(int i = 0; i < size; i++)
{
if(v(i) < minc)