]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add failing test on intel 15 607/head
authorTimo Heister <timo.heister@gmail.com>
Wed, 25 Feb 2015 22:33:07 +0000 (17:33 -0500)
committerTimo Heister <timo.heister@gmail.com>
Wed, 25 Feb 2015 22:34:58 +0000 (17:34 -0500)
This adds a simplified version of vector-vector.cc that fails in release mode on Intel.

tests/lac/intel-15-bug.cc [new file with mode: 0644]
tests/lac/intel-15-bug.output [new file with mode: 0644]

diff --git a/tests/lac/intel-15-bug.cc b/tests/lac/intel-15-bug.cc
new file mode 100644 (file)
index 0000000..ba9c68f
--- /dev/null
@@ -0,0 +1,90 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 1998 - 2015 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+// intel 15.0 compiler bug. This is a simplification of tests/lac/vector-vector
+// It only triggers when using TBB (it will use two threads) and only with
+// SIMD for long double.
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/lac/vector.h>
+#include <cmath>
+#include <fstream>
+#include <iomanip>
+#include <deal.II/base/parallel.h>
+
+typedef int size_type;
+
+template <typename Number>
+  struct Vectorization_add_v
+{
+  Number *val;
+  Number *v_val;
+
+  void operator() (const tbb::blocked_range<size_type> &range) const
+  {
+      DEAL_II_OPENMP_SIMD_PRAGMA
+        for (size_type i=range.begin(); i<range.end(); ++i)
+          val[i] += v_val[i];
+  }
+};
+
+const unsigned int N=3;
+
+void check()
+{
+  std::vector<long double> d1(N), d2(N);
+  for (unsigned int i=0; i<N; ++i)
+    {
+      d1[i] = 1.0;
+      d2[i] = i;
+    };
+
+  Vectorization_add_v<long double> vector_add;
+  vector_add.val = &d1[0];
+  vector_add.v_val = &d2[0];
+  
+  if (1)
+    {
+      //fails:
+      tbb::parallel_for (tbb::blocked_range<size_type> (0,
+                                                       N,
+                                                       2),
+                        vector_add,
+                        tbb::auto_partitioner());
+    }
+  else
+    {
+      // works:
+      vector_add(tbb::blocked_range<int> (0,1,2));
+      vector_add(tbb::blocked_range<int> (1,3,2));
+    }
+
+  deallog << "result: " << d1[N-1] << " should be 3" << std::endl;
+}
+
+int main()
+{
+  std::ofstream logfile("output");
+  deallog << std::fixed;
+  deallog << std::setprecision(2);
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+
+  check();
+}
+
+
diff --git a/tests/lac/intel-15-bug.output b/tests/lac/intel-15-bug.output
new file mode 100644 (file)
index 0000000..d92c1c6
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL::result: 3.00 should be 3

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.