]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix a bug where add_and_dot was giving the wrong results if the vector size was 7785/head
authorBruno Turcksin <bruno.turcksin@gmail.com>
Wed, 6 Mar 2019 22:48:18 +0000 (22:48 +0000)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Fri, 8 Mar 2019 13:55:07 +0000 (13:55 +0000)
greater than 4096

doc/news/changes/minor/20190306BrunoTurcksin [new file with mode: 0644]
include/deal.II/lac/cuda_kernels.templates.h
tests/cuda/cuda_vector_05.cu [new file with mode: 0644]
tests/cuda/cuda_vector_05.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20190306BrunoTurcksin b/doc/news/changes/minor/20190306BrunoTurcksin
new file mode 100644 (file)
index 0000000..5bbf2b5
--- /dev/null
@@ -0,0 +1,4 @@
+Fixed: LinearAlgebra::CUDAWrappers::Vector::add_and_dot was giving a wrong result
+if the size of vector was greater than 4096.
+<br>
+(Bruno Turcksin, 2018/03/06)
index cadeb4d11a39eda8708dc1097569def637bb4f27..d42e469329e188e369a658baec4e135b4042d8e2 100644 (file)
@@ -517,7 +517,7 @@ namespace LinearAlgebra
         else
           res_buf[local_idx] = 0.;
 
-        for (unsigned int i = 1; i < block_size; ++i)
+        for (unsigned int i = 1; i < chunk_size; ++i)
           {
             const unsigned int idx = global_idx + i * block_size;
             if (idx < N)
diff --git a/tests/cuda/cuda_vector_05.cu b/tests/cuda/cuda_vector_05.cu
new file mode 100644 (file)
index 0000000..d11c598
--- /dev/null
@@ -0,0 +1,71 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 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.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/lac/cuda_vector.h>
+#include <deal.II/lac/read_write_vector.h>
+
+#include <fstream>
+#include <iostream>
+#include <vector>
+
+#include "../tests.h"
+
+// There was a bug where add_and_dot would give the wrong result if the size of
+// the vectors was greater than BLOCK_SIZE*CHUNK_SIZE.
+
+void
+test()
+{
+  const unsigned int                          size = 4100;
+  LinearAlgebra::CUDAWrappers::Vector<double> a(size);
+  LinearAlgebra::CUDAWrappers::Vector<double> b(size);
+  Vector<double>                              a_host(size);
+  Vector<double>                              b_host(size);
+
+  LinearAlgebra::ReadWriteVector<double> read_write_1(size);
+  LinearAlgebra::ReadWriteVector<double> read_write_2(size);
+
+  for (unsigned int i = 0; i < size; ++i)
+    {
+      read_write_1[i] = i;
+      read_write_2[i] = 5. + i;
+      a_host[i]       = i;
+      b_host[i]       = 5. + i;
+    }
+
+  a.import(read_write_1, VectorOperation::insert);
+  b.import(read_write_2, VectorOperation::insert);
+  AssertThrow(a.add_and_dot(2., a, b) == a_host.add_and_dot(2., a_host, b_host),
+              ExcMessage("Problem in add_and_dot"));
+}
+
+
+int
+main(int argc, char **argv)
+{
+  initlog();
+  deallog.depth_console(0);
+
+  init_cuda();
+
+  test();
+
+  deallog << "OK" << std::endl;
+
+  return 0;
+}
diff --git a/tests/cuda/cuda_vector_05.output b/tests/cuda/cuda_vector_05.output
new file mode 100644 (file)
index 0000000..0fd8fc1
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL::OK

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.