aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/mixingtypes.cpp
diff options
context:
space:
mode:
authorGravatar Christoph Hertzberg <chtz@informatik.uni-bremen.de>2018-07-26 00:01:24 +0200
committerGravatar Christoph Hertzberg <chtz@informatik.uni-bremen.de>2018-07-26 00:01:24 +0200
commit397b0547e1e3151baa64c9677bef5882ad24d1ea (patch)
treefbc874e9db7f1ec85b5541a4199521f504755095 /test/mixingtypes.cpp
parent5e79402b4a742ef33574a568689c70be4f3d8549 (diff)
DIsable static assertions only when necessary and disable double-promotion warnings in that case as well
Diffstat (limited to 'test/mixingtypes.cpp')
-rw-r--r--test/mixingtypes.cpp58
1 files changed, 43 insertions, 15 deletions
diff --git a/test/mixingtypes.cpp b/test/mixingtypes.cpp
index 38f062f1e..aad63ec2b 100644
--- a/test/mixingtypes.cpp
+++ b/test/mixingtypes.cpp
@@ -8,13 +8,27 @@
// 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/.
-// work around "uninitialized" warnings and give that option some testing
-#define EIGEN_INITIALIZE_MATRICES_BY_ZERO
+#if defined(EIGEN_TEST_PART_7)
#ifndef EIGEN_NO_STATIC_ASSERT
#define EIGEN_NO_STATIC_ASSERT // turn static asserts into runtime asserts in order to check them
#endif
+// ignore double-promotion diagnostic for clang and gcc, if we check for static assertion anyway:
+// TODO do the same for MSVC?
+#if defined(__clang__)
+# if (__clang_major__ * 100 + __clang_minor__) >= 308
+# pragma clang diagnostic ignored "-Wdouble-promotion"
+# endif
+#elif defined(__GNUC__)
+ // TODO is there a minimal GCC version for this? At least g++-4.7 seems to be fine with this.
+# pragma GCC diagnostic ignored "-Wdouble-promotion"
+#endif
+
+#endif
+
+
+
#if defined(EIGEN_TEST_PART_1) || defined(EIGEN_TEST_PART_2) || defined(EIGEN_TEST_PART_3)
#ifndef EIGEN_DONT_VECTORIZE
@@ -35,6 +49,28 @@ using namespace std;
VERIFY_IS_APPROX(XPR,REF); \
VERIFY( g_called && #XPR" not properly optimized");
+template<int SizeAtCompileType>
+void raise_assertion(Index size = SizeAtCompileType)
+{
+ // VERIFY_RAISES_ASSERT(mf+md); // does not even compile
+ Matrix<float, SizeAtCompileType, 1> vf; vf.setRandom(size);
+ Matrix<double, SizeAtCompileType, 1> vd; vd.setRandom(size);
+ VERIFY_RAISES_ASSERT(vf=vd);
+ VERIFY_RAISES_ASSERT(vf+=vd);
+ VERIFY_RAISES_ASSERT(vf-=vd);
+ VERIFY_RAISES_ASSERT(vd=vf);
+ VERIFY_RAISES_ASSERT(vd+=vf);
+ VERIFY_RAISES_ASSERT(vd-=vf);
+
+ // vd.asDiagonal() * mf; // does not even compile
+ // vcd.asDiagonal() * mf; // does not even compile
+
+#if 0 // we get other compilation errors here than just static asserts
+ VERIFY_RAISES_ASSERT(vd.dot(vf));
+#endif
+}
+
+
template<int SizeAtCompileType> void mixingtypes(int size = SizeAtCompileType)
{
typedef std::complex<float> CF;
@@ -73,13 +109,6 @@ template<int SizeAtCompileType> void mixingtypes(int size = SizeAtCompileType)
while(std::abs(scf)<epsf) scf = internal::random<CF>();
while(std::abs(scd)<epsd) scd = internal::random<CD>();
-// VERIFY_RAISES_ASSERT(mf+md); // does not even compile
-
-#ifdef EIGEN_DONT_VECTORIZE
- VERIFY_RAISES_ASSERT(vf=vd);
- VERIFY_RAISES_ASSERT(vf+=vd);
-#endif
-
// check scalar products
VERIFY_MIX_SCALAR(vcf * sf , vcf * complex<float>(sf));
VERIFY_MIX_SCALAR(sd * vcd , complex<double>(sd) * vcd);
@@ -119,9 +148,6 @@ template<int SizeAtCompileType> void mixingtypes(int size = SizeAtCompileType)
// check dot product
vf.dot(vf);
-#if 0 // we get other compilation errors here than just static asserts
- VERIFY_RAISES_ASSERT(vd.dot(vf));
-#endif
VERIFY_IS_APPROX(vcf.dot(vf), vcf.dot(vf.template cast<complex<float> >()));
// check diagonal product
@@ -130,9 +156,6 @@ template<int SizeAtCompileType> void mixingtypes(int size = SizeAtCompileType)
VERIFY_IS_APPROX(mcf * vf.asDiagonal(), mcf * vf.template cast<complex<float> >().asDiagonal());
VERIFY_IS_APPROX(md * vcd.asDiagonal(), md.template cast<complex<double> >() * vcd.asDiagonal());
-// vd.asDiagonal() * mf; // does not even compile
-// vcd.asDiagonal() * mf; // does not even compile
-
// check inner product
VERIFY_IS_APPROX((vf.transpose() * vcf).value(), (vf.template cast<complex<float> >().transpose() * vcf).value());
@@ -296,5 +319,10 @@ EIGEN_DECLARE_TEST(mixingtypes)
CALL_SUBTEST_4(mixingtypes<3>());
CALL_SUBTEST_5(mixingtypes<4>());
CALL_SUBTEST_6(mixingtypes<Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE)));
+ CALL_SUBTEST_7(raise_assertion<Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE)));
}
+ CALL_SUBTEST_7(raise_assertion<0>());
+ CALL_SUBTEST_7(raise_assertion<3>());
+ CALL_SUBTEST_7(raise_assertion<4>());
+ CALL_SUBTEST_7(raise_assertion<Dynamic>(0));
}