aboutsummaryrefslogtreecommitdiff
path: root/test/nullary.cpp
diff options
context:
space:
mode:
Diffstat (limited to 'test/nullary.cpp')
-rw-r--r--test/nullary.cpp157
1 files changed, 97 insertions, 60 deletions
diff --git a/test/nullary.cpp b/test/nullary.cpp
index acd55506e..9b25ea4f3 100644
--- a/test/nullary.cpp
+++ b/test/nullary.cpp
@@ -70,7 +70,7 @@ void testVectorType(const VectorType& base)
Scalar high = internal::random<Scalar>(-500,500);
Scalar low = (size == 1 ? high : internal::random<Scalar>(-500,500));
- if (low>high) std::swap(low,high);
+ if (numext::real(low)>numext::real(high)) std::swap(low,high);
// check low==high
if(internal::random<float>(0.f,1.f)<0.05f)
@@ -79,7 +79,7 @@ void testVectorType(const VectorType& base)
else if(size>2 && std::numeric_limits<RealScalar>::max_exponent10>0 && internal::random<float>(0.f,1.f)<0.1f)
low = -internal::random<Scalar>(1,2) * RealScalar(std::pow(RealScalar(10),std::numeric_limits<RealScalar>::max_exponent10/2));
- const Scalar step = ((size == 1) ? 1 : (high-low)/(size-1));
+ const Scalar step = ((size == 1) ? 1 : (high-low)/RealScalar(size-1));
// check whether the result yields what we expect it to do
VectorType m(base);
@@ -89,21 +89,22 @@ void testVectorType(const VectorType& base)
{
VectorType n(size);
for (int i=0; i<size; ++i)
- n(i) = low+i*step;
+ n(i) = low+RealScalar(i)*step;
VERIFY_IS_APPROX(m,n);
CALL_SUBTEST( check_extremity_accuracy(m, low, high) );
}
- if((!NumTraits<Scalar>::IsInteger) || ((high-low)>=size && (Index(high-low)%(size-1))==0) || (Index(high-low+1)<size && (size%Index(high-low+1))==0))
+ RealScalar range_length = numext::real(high-low);
+ if((!NumTraits<Scalar>::IsInteger) || (range_length>=size && (Index(range_length)%(size-1))==0) || (Index(range_length+1)<size && (size%Index(range_length+1))==0))
{
VectorType n(size);
- if((!NumTraits<Scalar>::IsInteger) || (high-low>=size))
+ if((!NumTraits<Scalar>::IsInteger) || (range_length>=size))
for (int i=0; i<size; ++i)
- n(i) = size==1 ? low : (low + ((high-low)*Scalar(i))/(size-1));
+ n(i) = size==1 ? low : (low + ((high-low)*Scalar(i))/RealScalar(size-1));
else
for (int i=0; i<size; ++i)
- n(i) = size==1 ? low : low + Scalar((double(high-low+1)*double(i))/double(size));
+ n(i) = size==1 ? low : low + Scalar((double(range_length+1)*double(i))/double(size));
VERIFY_IS_APPROX(m,n);
// random access version
@@ -116,12 +117,12 @@ void testVectorType(const VectorType& base)
CALL_SUBTEST( check_extremity_accuracy(m, low, high) );
}
- VERIFY( m(m.size()-1) <= high );
- VERIFY( (m.array() <= high).all() );
- VERIFY( (m.array() >= low).all() );
+ VERIFY( numext::real(m(m.size()-1)) <= numext::real(high) );
+ VERIFY( (m.array().real() <= numext::real(high)).all() );
+ VERIFY( (m.array().real() >= numext::real(low)).all() );
- VERIFY( m(m.size()-1) >= low );
+ VERIFY( numext::real(m(m.size()-1)) >= numext::real(low) );
if(size>=1)
{
VERIFY( internal::isApprox(m(0),low) );
@@ -135,7 +136,7 @@ void testVectorType(const VectorType& base)
col_vector.setLinSpaced(size,low,high);
// when using the extended precision (e.g., FPU) the relative error might exceed 1 bit
// when computing the squared sum in isApprox, thus the 2x factor.
- VERIFY( row_vector.isApprox(col_vector.transpose(), Scalar(2)*NumTraits<Scalar>::epsilon()));
+ VERIFY( row_vector.isApprox(col_vector.transpose(), RealScalar(2)*NumTraits<Scalar>::epsilon()));
Matrix<Scalar,Dynamic,1> size_changer(size+50);
size_changer.setLinSpaced(size,low,high);
@@ -157,18 +158,18 @@ void testVectorType(const VectorType& base)
{
Index n0 = VectorType::SizeAtCompileTime==Dynamic ? 0 : VectorType::SizeAtCompileTime;
low = internal::random<Scalar>();
- m = VectorType::LinSpaced(n0,low,low-1);
+ m = VectorType::LinSpaced(n0,low,low-RealScalar(1));
VERIFY(m.size()==n0);
if(VectorType::SizeAtCompileTime==Dynamic)
{
VERIFY_IS_EQUAL(VectorType::LinSpaced(n0,0,Scalar(n0-1)).sum(),Scalar(0));
- VERIFY_IS_EQUAL(VectorType::LinSpaced(n0,low,low-1).sum(),Scalar(0));
+ VERIFY_IS_EQUAL(VectorType::LinSpaced(n0,low,low-RealScalar(1)).sum(),Scalar(0));
}
m.setLinSpaced(n0,0,Scalar(n0-1));
VERIFY(m.size()==n0);
- m.setLinSpaced(n0,low,low-1);
+ m.setLinSpaced(n0,low,low-RealScalar(1));
VERIFY(m.size()==n0);
// empty range only:
@@ -178,19 +179,37 @@ void testVectorType(const VectorType& base)
if(NumTraits<Scalar>::IsInteger)
{
- VERIFY_IS_APPROX( VectorType::LinSpaced(size,low,Scalar(low+size-1)), VectorType::LinSpaced(size,Scalar(low+size-1),low).reverse() );
+ VERIFY_IS_APPROX( VectorType::LinSpaced(size,low,low+Scalar(size-1)), VectorType::LinSpaced(size,low+Scalar(size-1),low).reverse() );
if(VectorType::SizeAtCompileTime==Dynamic)
{
// Check negative multiplicator path:
for(Index k=1; k<5; ++k)
- VERIFY_IS_APPROX( VectorType::LinSpaced(size,low,Scalar(low+(size-1)*k)), VectorType::LinSpaced(size,Scalar(low+(size-1)*k),low).reverse() );
+ VERIFY_IS_APPROX( VectorType::LinSpaced(size,low,low+Scalar((size-1)*k)), VectorType::LinSpaced(size,low+Scalar((size-1)*k),low).reverse() );
// Check negative divisor path:
for(Index k=1; k<5; ++k)
- VERIFY_IS_APPROX( VectorType::LinSpaced(size*k,low,Scalar(low+size-1)), VectorType::LinSpaced(size*k,Scalar(low+size-1),low).reverse() );
+ VERIFY_IS_APPROX( VectorType::LinSpaced(size*k,low,low+Scalar(size-1)), VectorType::LinSpaced(size*k,low+Scalar(size-1),low).reverse() );
}
}
}
+
+ // test setUnit()
+ if(m.size()>0)
+ {
+ for(Index k=0; k<10; ++k)
+ {
+ Index i = internal::random<Index>(0,m.size()-1);
+ m.setUnit(i);
+ VERIFY_IS_APPROX( m, VectorType::Unit(m.size(), i) );
+ }
+ if(VectorType::SizeAtCompileTime==Dynamic)
+ {
+ Index i = internal::random<Index>(0,2*m.size()-1);
+ m.setUnit(2*m.size(),i);
+ VERIFY_IS_APPROX( m, VectorType::Unit(m.size(),i) );
+ }
+ }
+
}
template<typename MatrixType>
@@ -221,45 +240,36 @@ void testMatrixType(const MatrixType& m)
VERIFY_IS_APPROX( A(i,j), s1 );
}
-void test_nullary()
+template<int>
+void bug79()
{
- CALL_SUBTEST_1( testMatrixType(Matrix2d()) );
- CALL_SUBTEST_2( testMatrixType(MatrixXcf(internal::random<int>(1,300),internal::random<int>(1,300))) );
- CALL_SUBTEST_3( testMatrixType(MatrixXf(internal::random<int>(1,300),internal::random<int>(1,300))) );
-
- for(int i = 0; i < g_repeat*10; i++) {
- CALL_SUBTEST_4( testVectorType(VectorXd(internal::random<int>(1,30000))) );
- CALL_SUBTEST_5( testVectorType(Vector4d()) ); // regression test for bug 232
- CALL_SUBTEST_6( testVectorType(Vector3d()) );
- CALL_SUBTEST_7( testVectorType(VectorXf(internal::random<int>(1,30000))) );
- CALL_SUBTEST_8( testVectorType(Vector3f()) );
- CALL_SUBTEST_8( testVectorType(Vector4f()) );
- CALL_SUBTEST_8( testVectorType(Matrix<float,8,1>()) );
- CALL_SUBTEST_8( testVectorType(Matrix<float,1,1>()) );
-
- CALL_SUBTEST_9( testVectorType(VectorXi(internal::random<int>(1,10))) );
- CALL_SUBTEST_9( testVectorType(VectorXi(internal::random<int>(9,300))) );
- CALL_SUBTEST_9( testVectorType(Matrix<int,1,1>()) );
- }
-
-#ifdef EIGEN_TEST_PART_6
// Assignment of a RowVectorXd to a MatrixXd (regression test for bug #79).
VERIFY( (MatrixXd(RowVectorXd::LinSpaced(3, 0, 1)) - RowVector3d(0, 0.5, 1)).norm() < std::numeric_limits<double>::epsilon() );
-#endif
+}
-#ifdef EIGEN_TEST_PART_9
+template<int>
+void bug1630()
+{
+ Array4d x4 = Array4d::LinSpaced(0.0, 1.0);
+ Array3d x3(Array4d::LinSpaced(0.0, 1.0).head(3));
+ VERIFY_IS_APPROX(x4.head(3), x3);
+}
+
+template<int>
+void nullary_overflow()
+{
// Check possible overflow issue
- {
- int n = 60000;
- ArrayXi a1(n), a2(n);
- a1.setLinSpaced(n, 0, n-1);
- for(int i=0; i<n; ++i)
- a2(i) = i;
- VERIFY_IS_APPROX(a1,a2);
- }
-#endif
+ int n = 60000;
+ ArrayXi a1(n), a2(n);
+ a1.setLinSpaced(n, 0, n-1);
+ for(int i=0; i<n; ++i)
+ a2(i) = i;
+ VERIFY_IS_APPROX(a1,a2);
+}
-#ifdef EIGEN_TEST_PART_10
+template<int>
+void nullary_internal_logic()
+{
// check some internal logic
VERIFY(( internal::has_nullary_operator<internal::scalar_constant_op<double> >::value ));
VERIFY(( !internal::has_unary_operator<internal::scalar_constant_op<double> >::value ));
@@ -271,10 +281,10 @@ void test_nullary()
VERIFY(( internal::has_binary_operator<internal::scalar_identity_op<double> >::value ));
VERIFY(( !internal::functor_has_linear_access<internal::scalar_identity_op<double> >::ret ));
- VERIFY(( !internal::has_nullary_operator<internal::linspaced_op<float,float> >::value ));
- VERIFY(( internal::has_unary_operator<internal::linspaced_op<float,float> >::value ));
- VERIFY(( !internal::has_binary_operator<internal::linspaced_op<float,float> >::value ));
- VERIFY(( internal::functor_has_linear_access<internal::linspaced_op<float,float> >::ret ));
+ VERIFY(( !internal::has_nullary_operator<internal::linspaced_op<float> >::value ));
+ VERIFY(( internal::has_unary_operator<internal::linspaced_op<float> >::value ));
+ VERIFY(( !internal::has_binary_operator<internal::linspaced_op<float> >::value ));
+ VERIFY(( internal::functor_has_linear_access<internal::linspaced_op<float> >::ret ));
// Regression unit test for a weird MSVC bug.
// Search "nullary_wrapper_workaround_msvc" in CoreEvaluators.h for the details.
@@ -295,10 +305,37 @@ void test_nullary()
VERIFY(( !internal::has_binary_operator<internal::scalar_constant_op<float> >::value ));
VERIFY(( internal::functor_has_linear_access<internal::scalar_constant_op<float> >::ret ));
- VERIFY(( !internal::has_nullary_operator<internal::linspaced_op<int,int> >::value ));
- VERIFY(( internal::has_unary_operator<internal::linspaced_op<int,int> >::value ));
- VERIFY(( !internal::has_binary_operator<internal::linspaced_op<int,int> >::value ));
- VERIFY(( internal::functor_has_linear_access<internal::linspaced_op<int,int> >::ret ));
+ VERIFY(( !internal::has_nullary_operator<internal::linspaced_op<int> >::value ));
+ VERIFY(( internal::has_unary_operator<internal::linspaced_op<int> >::value ));
+ VERIFY(( !internal::has_binary_operator<internal::linspaced_op<int> >::value ));
+ VERIFY(( internal::functor_has_linear_access<internal::linspaced_op<int> >::ret ));
}
-#endif
+}
+
+EIGEN_DECLARE_TEST(nullary)
+{
+ CALL_SUBTEST_1( testMatrixType(Matrix2d()) );
+ CALL_SUBTEST_2( testMatrixType(MatrixXcf(internal::random<int>(1,300),internal::random<int>(1,300))) );
+ CALL_SUBTEST_3( testMatrixType(MatrixXf(internal::random<int>(1,300),internal::random<int>(1,300))) );
+
+ for(int i = 0; i < g_repeat*10; i++) {
+ CALL_SUBTEST_3( testVectorType(VectorXcd(internal::random<int>(1,30000))) );
+ CALL_SUBTEST_4( testVectorType(VectorXd(internal::random<int>(1,30000))) );
+ CALL_SUBTEST_5( testVectorType(Vector4d()) ); // regression test for bug 232
+ CALL_SUBTEST_6( testVectorType(Vector3d()) );
+ CALL_SUBTEST_7( testVectorType(VectorXf(internal::random<int>(1,30000))) );
+ CALL_SUBTEST_8( testVectorType(Vector3f()) );
+ CALL_SUBTEST_8( testVectorType(Vector4f()) );
+ CALL_SUBTEST_8( testVectorType(Matrix<float,8,1>()) );
+ CALL_SUBTEST_8( testVectorType(Matrix<float,1,1>()) );
+
+ CALL_SUBTEST_9( testVectorType(VectorXi(internal::random<int>(1,10))) );
+ CALL_SUBTEST_9( testVectorType(VectorXi(internal::random<int>(9,300))) );
+ CALL_SUBTEST_9( testVectorType(Matrix<int,1,1>()) );
+ }
+
+ CALL_SUBTEST_6( bug79<0>() );
+ CALL_SUBTEST_6( bug1630<0>() );
+ CALL_SUBTEST_9( nullary_overflow<0>() );
+ CALL_SUBTEST_10( nullary_internal_logic<0>() );
}