aboutsummaryrefslogtreecommitdiffhomepage
path: root/test
diff options
context:
space:
mode:
Diffstat (limited to 'test')
-rw-r--r--test/CMakeLists.txt3
-rw-r--r--test/block.cpp10
-rw-r--r--test/cholesky.cpp24
-rw-r--r--test/mapstaticmethods.cpp6
-rw-r--r--test/packetmath.cpp15
-rw-r--r--test/product_notemporary.cpp5
-rw-r--r--test/ref.cpp5
-rw-r--r--test/sparse.h4
-rw-r--r--test/sparse_basic.cpp8
-rw-r--r--test/triangular.cpp7
10 files changed, 59 insertions, 28 deletions
diff --git a/test/CMakeLists.txt b/test/CMakeLists.txt
index 8321328cd..154e62424 100644
--- a/test/CMakeLists.txt
+++ b/test/CMakeLists.txt
@@ -82,7 +82,7 @@ endif()
find_package(Pastix)
find_package(Scotch)
-find_package(Metis)
+find_package(Metis 5.0 REQUIRED)
if(PASTIX_FOUND)
add_definitions("-DEIGEN_PASTIX_SUPPORT")
include_directories(${PASTIX_INCLUDES})
@@ -304,6 +304,7 @@ ei_add_property(EIGEN_TESTING_SUMMARY "CXX_FLAGS: ${CMAKE_CXX_FLAGS}\n")
ei_add_property(EIGEN_TESTING_SUMMARY "Sparse lib flags: ${SPARSE_LIBS}\n")
option(EIGEN_TEST_EIGEN2 "Run whole Eigen2 test suite against EIGEN2_SUPPORT" OFF)
+mark_as_advanced(EIGEN_TEST_EIGEN2)
if(EIGEN_TEST_EIGEN2)
message(WARNING "The Eigen2 test suite has been removed")
endif()
diff --git a/test/block.cpp b/test/block.cpp
index 9ed5d7bc5..269acd28e 100644
--- a/test/block.cpp
+++ b/test/block.cpp
@@ -141,11 +141,11 @@ template<typename MatrixType> void block(const MatrixType& m)
VERIFY_IS_EQUAL( (m1.transpose().block(c1,r1,c2-c1+1,r2-r1+1).col(0)) , (m1.row(r1).segment(c1,c2-c1+1)).transpose() );
// expressions without direct access
- VERIFY_IS_EQUAL( ((m1+m2).block(r1,c1,rows-r1,cols-c1).block(r2-r1,c2-c1,rows-r2,cols-c2)) , ((m1+m2).block(r2,c2,rows-r2,cols-c2)) );
- VERIFY_IS_EQUAL( ((m1+m2).block(r1,c1,r2-r1+1,c2-c1+1).row(0)) , ((m1+m2).row(r1).segment(c1,c2-c1+1)) );
- VERIFY_IS_EQUAL( ((m1+m2).block(r1,c1,r2-r1+1,c2-c1+1).col(0)) , ((m1+m2).col(c1).segment(r1,r2-r1+1)) );
- VERIFY_IS_EQUAL( ((m1+m2).block(r1,c1,r2-r1+1,c2-c1+1).transpose().col(0)) , ((m1+m2).row(r1).segment(c1,c2-c1+1)).transpose() );
- VERIFY_IS_EQUAL( ((m1+m2).transpose().block(c1,r1,c2-c1+1,r2-r1+1).col(0)) , ((m1+m2).row(r1).segment(c1,c2-c1+1)).transpose() );
+ VERIFY_IS_APPROX( ((m1+m2).block(r1,c1,rows-r1,cols-c1).block(r2-r1,c2-c1,rows-r2,cols-c2)) , ((m1+m2).block(r2,c2,rows-r2,cols-c2)) );
+ VERIFY_IS_APPROX( ((m1+m2).block(r1,c1,r2-r1+1,c2-c1+1).row(0)) , ((m1+m2).row(r1).segment(c1,c2-c1+1)) );
+ VERIFY_IS_APPROX( ((m1+m2).block(r1,c1,r2-r1+1,c2-c1+1).col(0)) , ((m1+m2).col(c1).segment(r1,r2-r1+1)) );
+ VERIFY_IS_APPROX( ((m1+m2).block(r1,c1,r2-r1+1,c2-c1+1).transpose().col(0)) , ((m1+m2).row(r1).segment(c1,c2-c1+1)).transpose() );
+ VERIFY_IS_APPROX( ((m1+m2).transpose().block(c1,r1,c2-c1+1,r2-r1+1).col(0)) , ((m1+m2).row(r1).segment(c1,c2-c1+1)).transpose() );
// evaluation into plain matrices from expressions with direct access (stress MapBase)
DynamicMatrixType dm;
diff --git a/test/cholesky.cpp b/test/cholesky.cpp
index 569318f83..a883192ab 100644
--- a/test/cholesky.cpp
+++ b/test/cholesky.cpp
@@ -68,6 +68,7 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
Index cols = m.cols();
typedef typename MatrixType::Scalar Scalar;
+ typedef typename NumTraits<Scalar>::Real RealScalar;
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
@@ -180,7 +181,7 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
if(rows>=3)
{
SquareMatrixType A = symm;
- int c = internal::random<int>(0,rows-2);
+ Index c = internal::random<Index>(0,rows-2);
A.bottomRightCorner(c,c).setZero();
// Make sure a solution exists:
vecX.setRandom();
@@ -195,7 +196,7 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
// check non-full rank matrices
if(rows>=3)
{
- int r = internal::random<int>(1,rows-1);
+ Index r = internal::random<Index>(1,rows-1);
Matrix<Scalar,Dynamic,Dynamic> a = Matrix<Scalar,Dynamic,Dynamic>::Random(rows,r);
SquareMatrixType A = a * a.adjoint();
// Make sure a solution exists:
@@ -207,6 +208,25 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
vecX = ldltlo.solve(vecB);
VERIFY_IS_APPROX(A * vecX, vecB);
}
+
+ // check matrices with a wide spectrum
+ if(rows>=3)
+ {
+ RealScalar s = (std::min)(16,std::numeric_limits<RealScalar>::max_exponent10/8);
+ Matrix<Scalar,Dynamic,Dynamic> a = Matrix<Scalar,Dynamic,Dynamic>::Random(rows,rows);
+ Matrix<RealScalar,Dynamic,1> d = Matrix<RealScalar,Dynamic,1>::Random(rows);
+ for(Index k=0; k<rows; ++k)
+ d(k) = d(k)*std::pow(RealScalar(10),internal::random<RealScalar>(-s,s));
+ SquareMatrixType A = a * d.asDiagonal() * a.adjoint();
+ // Make sure a solution exists:
+ vecX.setRandom();
+ vecB = A * vecX;
+ vecX.setZero();
+ ldltlo.compute(A);
+ VERIFY_IS_APPROX(A, ldltlo.reconstructedMatrix());
+ vecX = ldltlo.solve(vecB);
+ VERIFY_IS_APPROX(A * vecX, vecB);
+ }
}
// update/downdate
diff --git a/test/mapstaticmethods.cpp b/test/mapstaticmethods.cpp
index 5b512bde4..06272d106 100644
--- a/test/mapstaticmethods.cpp
+++ b/test/mapstaticmethods.cpp
@@ -69,7 +69,8 @@ struct mapstaticmethods_impl<PlainObjectType, true, false>
{
static void run(const PlainObjectType& m)
{
- int rows = m.rows(), cols = m.cols();
+ typedef typename PlainObjectType::Index Index;
+ Index rows = m.rows(), cols = m.cols();
int i = internal::random<int>(2,5), j = internal::random<int>(2,5);
@@ -115,7 +116,8 @@ struct mapstaticmethods_impl<PlainObjectType, true, true>
{
static void run(const PlainObjectType& v)
{
- int size = v.size();
+ typedef typename PlainObjectType::Index Index;
+ Index size = v.size();
int i = internal::random<int>(2,5);
diff --git a/test/packetmath.cpp b/test/packetmath.cpp
index a51d31dbd..e5dc473c2 100644
--- a/test/packetmath.cpp
+++ b/test/packetmath.cpp
@@ -408,14 +408,17 @@ template<typename Scalar> void packetmath_scatter_gather() {
for (int i=0; i<PacketSize; ++i) {
data1[i] = internal::random<Scalar>()/RealScalar(PacketSize);
}
- EIGEN_ALIGN_DEFAULT Scalar buffer[PacketSize*11];
- memset(buffer, 0, 11*sizeof(Packet));
+
+ int stride = internal::random<int>(1,20);
+
+ EIGEN_ALIGN_DEFAULT Scalar buffer[PacketSize*20];
+ memset(buffer, 0, 20*sizeof(Packet));
Packet packet = internal::pload<Packet>(data1);
- internal::pscatter<Scalar, Packet>(buffer, packet, 11);
+ internal::pscatter<Scalar, Packet>(buffer, packet, stride);
- for (int i = 0; i < PacketSize*11; ++i) {
- if ((i%11) == 0) {
- VERIFY(isApproxAbs(buffer[i], data1[i/11], refvalue) && "pscatter");
+ for (int i = 0; i < PacketSize*20; ++i) {
+ if ((i%stride) == 0 && i<stride*PacketSize) {
+ VERIFY(isApproxAbs(buffer[i], data1[i/stride], refvalue) && "pscatter");
} else {
VERIFY(isApproxAbs(buffer[i], Scalar(0), refvalue) && "pscatter");
}
diff --git a/test/product_notemporary.cpp b/test/product_notemporary.cpp
index deb494c96..805cc8939 100644
--- a/test/product_notemporary.cpp
+++ b/test/product_notemporary.cpp
@@ -7,13 +7,12 @@
// 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/.
-static int nb_temporaries;
+static long int nb_temporaries;
-inline void on_temporary_creation(int size) {
+inline void on_temporary_creation(long int size) {
// here's a great place to set a breakpoint when debugging failures in this test!
if(size!=0) nb_temporaries++;
}
-
#define EIGEN_DENSE_STORAGE_CTOR_PLUGIN { on_temporary_creation(size); }
diff --git a/test/ref.cpp b/test/ref.cpp
index 19e81549c..d91e3b54c 100644
--- a/test/ref.cpp
+++ b/test/ref.cpp
@@ -12,13 +12,12 @@
#undef EIGEN_DEFAULT_TO_ROW_MAJOR
#endif
-static int nb_temporaries;
+static long int nb_temporaries;
-inline void on_temporary_creation(int) {
+inline void on_temporary_creation(long int) {
// here's a great place to set a breakpoint when debugging failures in this test!
nb_temporaries++;
}
-
#define EIGEN_DENSE_STORAGE_CTOR_PLUGIN { on_temporary_creation(size); }
diff --git a/test/sparse.h b/test/sparse.h
index e19a76316..81ab9e873 100644
--- a/test/sparse.h
+++ b/test/sparse.h
@@ -71,7 +71,7 @@ initSparse(double density,
//sparseMat.startVec(j);
for(Index i=0; i<sparseMat.innerSize(); i++)
{
- int ai(i), aj(j);
+ Index ai(i), aj(j);
if(IsRowMajor)
std::swap(ai,aj);
Scalar v = (internal::random<double>(0,1) < density) ? internal::random<Scalar>() : Scalar(0);
@@ -163,7 +163,7 @@ initSparse(double density,
{
sparseVec.reserve(int(refVec.size()*density));
sparseVec.setZero();
- for(Index i=0; i<refVec.size(); i++)
+ for(int i=0; i<refVec.size(); i++)
{
Scalar v = (internal::random<double>(0,1) < density) ? internal::random<Scalar>() : Scalar(0);
if (v!=Scalar(0))
diff --git a/test/sparse_basic.cpp b/test/sparse_basic.cpp
index 6c620f037..4c9b9111e 100644
--- a/test/sparse_basic.cpp
+++ b/test/sparse_basic.cpp
@@ -147,7 +147,7 @@ template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& re
DenseMatrix m1(rows,cols);
m1.setZero();
SparseMatrixType m2(rows,cols);
- VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? m2.innerSize() : std::max<int>(1,m2.innerSize()/8)));
+ VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? int(m2.innerSize()) : std::max<int>(1,int(m2.innerSize())/8)));
m2.reserve(r);
for (int k=0; k<rows*cols; ++k)
{
@@ -181,7 +181,7 @@ template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& re
VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.col(j0)+refMat2.col(j1));
SparseMatrixType m3(rows,rows);
- m3.reserve(VectorXi::Constant(rows,rows/2));
+ m3.reserve(VectorXi::Constant(rows,int(rows/2)));
for(Index j=0; j<rows; ++j)
for(Index k=0; k<j; ++k)
m3.insertByOuterInner(j,k) = k+1;
@@ -384,11 +384,11 @@ template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& re
{
typedef Triplet<Scalar,Index> TripletType;
std::vector<TripletType> triplets;
- int ntriplets = rows*cols;
+ Index ntriplets = rows*cols;
triplets.reserve(ntriplets);
DenseMatrix refMat(rows,cols);
refMat.setZero();
- for(int i=0;i<ntriplets;++i)
+ for(Index i=0;i<ntriplets;++i)
{
Index r = internal::random<Index>(0,rows-1);
Index c = internal::random<Index>(0,cols-1);
diff --git a/test/triangular.cpp b/test/triangular.cpp
index 54320390b..936c2aef3 100644
--- a/test/triangular.cpp
+++ b/test/triangular.cpp
@@ -113,6 +113,13 @@ template<typename MatrixType> void triangular_square(const MatrixType& m)
m3.setZero();
m3.template triangularView<Upper>().setOnes();
VERIFY_IS_APPROX(m2,m3);
+
+ m1.setRandom();
+ m3 = m1.template triangularView<Upper>();
+ Matrix<Scalar, MatrixType::ColsAtCompileTime, Dynamic> m5(cols, internal::random<int>(1,20)); m5.setRandom();
+ Matrix<Scalar, Dynamic, MatrixType::RowsAtCompileTime> m6(internal::random<int>(1,20), rows); m6.setRandom();
+ VERIFY_IS_APPROX(m1.template triangularView<Upper>() * m5, m3*m5);
+ VERIFY_IS_APPROX(m6*m1.template triangularView<Upper>(), m6*m3);
}