]> https://gitweb.dealii.org/ - dealii.git/commitdiff
extend internal::VectorOperations::parallel_for() to [start,end) 4015/head
authorDenis Davydov <davydden@gmail.com>
Wed, 1 Mar 2017 07:11:30 +0000 (08:11 +0100)
committerDenis Davydov <davydden@gmail.com>
Wed, 1 Mar 2017 09:47:33 +0000 (10:47 +0100)
include/deal.II/lac/vector_operations_internal.h
tests/lac/vector_operations_parallel_for_start_end.cc [new file with mode: 0644]
tests/lac/vector_operations_parallel_for_start_end.output [new file with mode: 0644]

index 57a2cf0601078ed6e29ae8abd2159a059519c961..2947443272ef7eee39cd3855d0b96ee0f8e37537 100644 (file)
@@ -116,11 +116,14 @@ namespace internal
     struct TBBForFunctor
     {
       TBBForFunctor(Functor &functor,
-                    const size_type vec_size)
+                    const size_type start,
+                    const size_type end)
         :
         functor(functor),
-        vec_size(vec_size)
+        start(start),
+        end(end)
       {
+        const size_type vec_size = end-start;
         // set chunk size for sub-tasks
         const unsigned int gs = internal::Vector::minimum_parallel_grain_size;
         n_chunks = std::min(static_cast<size_type>(4*MultithreadInfo::n_threads()),
@@ -140,13 +143,14 @@ namespace internal
 
       void operator() (const tbb::blocked_range<size_type> &range) const
       {
-        const size_type begin = range.begin()*chunk_size;
-        const size_type end = std::min(range.end()*chunk_size, vec_size);
-        functor(begin, end);
+        const size_type r_begin = start + range.begin()*chunk_size;
+        const size_type r_end = std::min(start + range.end()*chunk_size, end);
+        functor(r_begin, r_end);
       }
 
       Functor &functor;
-      const size_type vec_size;
+      const size_type start;
+      const size_type end;
       unsigned int n_chunks;
       size_type chunk_size;
     };
@@ -154,9 +158,11 @@ namespace internal
 
     template <typename Functor>
     void parallel_for(Functor &functor,
-                      size_type vec_size,
-                      std_cxx11::shared_ptr<parallel::internal::TBBPartitioner> &partitioner)
+                      size_type end,
+                      std_cxx11::shared_ptr<parallel::internal::TBBPartitioner> &partitioner,
+                      size_type start = 0)
     {
+      size_type vec_size = end-start;
 #ifdef DEAL_II_WITH_THREADS
       // only go to the parallel function in case there are at least 4 parallel
       // items, otherwise the overhead is too large
@@ -169,7 +175,7 @@ namespace internal
           std_cxx11::shared_ptr<tbb::affinity_partitioner> tbb_partitioner =
             partitioner->acquire_one_partitioner();
 
-          TBBForFunctor<Functor> generic_functor(functor, vec_size);
+          TBBForFunctor<Functor> generic_functor(functor, start, end);
           tbb::parallel_for (tbb::blocked_range<size_type> (0,
                                                             generic_functor.n_chunks,
                                                             1),
@@ -178,9 +184,9 @@ namespace internal
           partitioner->release_one_partitioner(tbb_partitioner);
         }
       else if (vec_size > 0)
-        functor(0,vec_size);
+        functor(start,end);
 #else
-      functor(0,vec_size);
+      functor(start,end);
       (void)partitioner;
 #endif
     }
diff --git a/tests/lac/vector_operations_parallel_for_start_end.cc b/tests/lac/vector_operations_parallel_for_start_end.cc
new file mode 100644 (file)
index 0000000..d7965db
--- /dev/null
@@ -0,0 +1,80 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2012 - 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.
+//
+// ---------------------------------------------------------------------
+
+
+// check that internal::VectorOperations::parallel_for works for start-end
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/index_set.h>
+#include <deal.II/lac/vector_operations_internal.h>
+#include <cmath>
+#include <fstream>
+#include <iomanip>
+
+
+
+
+template <typename Number>
+void check ()
+{
+  for (unsigned int test=0; test<5; ++test)
+    {
+      const unsigned int size = 17 + test*1101;
+
+      std_cxx11::shared_ptr< ::dealii::parallel::internal::TBBPartitioner> thread_loop_partitioner;
+      thread_loop_partitioner.reset(new ::dealii::parallel::internal::TBBPartitioner());
+
+      Number *val;
+      Utilities::System::posix_memalign ((void **)&val, 64, sizeof(Number)*size);
+
+      const Number s = 3.1415;
+      internal::VectorOperations::Vector_set<Number> setter(s, val);
+
+      // now break the size in chunks
+      const unsigned int n_chunks = 3;
+      const unsigned int chunk_size = size / n_chunks;
+      for (unsigned int i = 0; i <= n_chunks; ++i)
+        {
+          const unsigned int begin = i*chunk_size;
+          const unsigned int end   = std::min((i+1)*chunk_size, size);
+          internal::VectorOperations::parallel_for(setter, end,
+                                                   thread_loop_partitioner,
+                                                   begin);
+        }
+
+      // check values:
+      for (unsigned int i = 0; i < size; ++i)
+        AssertThrow(val[i] == s,
+                    ExcInternalError());
+
+      free(val);
+    }
+}
+
+
+int main()
+{
+  std::ofstream logfile("output");
+  deallog << std::fixed;
+  deallog << std::setprecision(2);
+  deallog.attach(logfile);
+  deallog.threshold_double(1.e-10);
+
+  check<float>();
+  check<double>();
+  check<long double>();
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/lac/vector_operations_parallel_for_start_end.output b/tests/lac/vector_operations_parallel_for_start_end.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.