# include <deal.II/base/config.h>
# include <deal.II/base/graph_coloring.h>
+# include <deal.II/base/iterator_range.h>
# include <deal.II/base/multithread_info.h>
# include <deal.II/base/parallel.h>
# include <deal.II/base/template_constraints.h>
# endif
# include <functional>
+# include <iterator>
# include <memory>
# include <utility>
# include <vector>
}
+
+ /**
+ * Same as the function above, but for iterator ranges and C-style arrays.
+ * A class that fulfills the requirements of an iterator range defines the
+ * functions `IteratorRangeType::begin()` and `IteratorRangeType::end()`,
+ * both of which return iterators to elements that form the bounds of the
+ * range.
+ */
+ template <typename Worker,
+ typename Copier,
+ typename IteratorRangeType,
+ typename ScratchData,
+ typename CopyData,
+ typename = typename std::enable_if<
+ has_begin_and_end<IteratorRangeType>::value>::type>
+ void
+ run(IteratorRangeType iterator_range,
+ Worker worker,
+ Copier copier,
+ const ScratchData &sample_scratch_data,
+ const CopyData & sample_copy_data,
+ const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
+ const unsigned int chunk_size = 8)
+ {
+ // Call the function above
+ run(iterator_range.begin(),
+ iterator_range.end(),
+ worker,
+ copier,
+ sample_scratch_data,
+ sample_copy_data,
+ queue_length,
+ chunk_size);
+ }
+
+
+
+ /**
+ * Same as the function above, but for deal.II's IteratorRange.
+ */
+ template <typename Worker,
+ typename Copier,
+ typename Iterator,
+ typename ScratchData,
+ typename CopyData>
+ void
+ run(const IteratorRange<Iterator> &iterator_range,
+ Worker worker,
+ Copier copier,
+ const ScratchData & sample_scratch_data,
+ const CopyData & sample_copy_data,
+ const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
+ const unsigned int chunk_size = 8)
+ {
+ // Call the function above
+ run(iterator_range.begin(),
+ iterator_range.end(),
+ worker,
+ copier,
+ sample_scratch_data,
+ sample_copy_data,
+ queue_length,
+ chunk_size);
+ }
+
+
// Implementation 3:
template <typename Worker,
typename Copier,
chunk_size);
}
+
+ template <typename MainClass,
+ typename Iterator,
+ typename ScratchData,
+ typename CopyData>
+ void
+ run(const typename IteratorRange<Iterator>::IteratorOverIterators &begin,
+ const typename IteratorRange<
+ typename identity<Iterator>::type>::IteratorOverIterators &end,
+ MainClass & main_object,
+ void (MainClass::*worker)(const Iterator &, ScratchData &, CopyData &),
+ void (MainClass::*copier)(const CopyData &),
+ const ScratchData &sample_scratch_data,
+ const CopyData & sample_copy_data,
+ const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
+ const unsigned int chunk_size = 8)
+ {
+ // forward to the other function
+ run(begin,
+ end,
+ std::bind(worker,
+ std::ref(main_object),
+ std::placeholders::_1,
+ std::placeholders::_2,
+ std::placeholders::_3),
+ std::bind(copier, std::ref(main_object), std::placeholders::_1),
+ sample_scratch_data,
+ sample_copy_data,
+ queue_length,
+ chunk_size);
+ }
+
+
+
+ /**
+ * Same as the function above, but for iterator ranges and C-style arrays.
+ * A class that fulfills the requirements of an iterator range defines the
+ * functions `IteratorRangeType::begin()` and `IteratorRangeType::end()`,
+ * both of which return iterators to elements that form the bounds of the
+ * range.
+ */
+ template <typename MainClass,
+ typename IteratorRangeType,
+ typename ScratchData,
+ typename CopyData,
+ typename = typename std::enable_if<
+ has_begin_and_end<IteratorRangeType>::value>::type>
+ void
+ run(IteratorRangeType iterator_range,
+ MainClass & main_object,
+ void (MainClass::*worker)(
+ const typename identity<IteratorRangeType>::type::iterator &,
+ ScratchData &,
+ CopyData &),
+ void (MainClass::*copier)(const CopyData &),
+ const ScratchData &sample_scratch_data,
+ const CopyData & sample_copy_data,
+ const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
+ const unsigned int chunk_size = 8)
+ {
+ // Call the function above
+ run(std::begin(iterator_range),
+ std::end(iterator_range),
+ main_object,
+ worker,
+ copier,
+ sample_scratch_data,
+ sample_copy_data,
+ queue_length,
+ chunk_size);
+ }
+
+
+
+ /**
+ * Same as the function above, but for deal.II's IteratorRange.
+ */
+ template <typename MainClass,
+ typename Iterator,
+ typename ScratchData,
+ typename CopyData>
+ void
+ run(IteratorRange<Iterator> iterator_range,
+ MainClass & main_object,
+ void (MainClass::*worker)(const Iterator &, ScratchData &, CopyData &),
+ void (MainClass::*copier)(const CopyData &),
+ const ScratchData &sample_scratch_data,
+ const CopyData & sample_copy_data,
+ const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
+ const unsigned int chunk_size = 8)
+ {
+ // Call the function above
+ run(std::begin(iterator_range),
+ std::end(iterator_range),
+ main_object,
+ worker,
+ copier,
+ sample_scratch_data,
+ sample_copy_data,
+ queue_length,
+ chunk_size);
+ }
+
} // namespace WorkStream
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2008 - 2017 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.
+//
+// ---------------------------------------------------------------------
+
+
+// test functions in namespace WorkStream
+// This test is based off base/work_stream_01.cc, but for iterator ranges
+
+#include <deal.II/base/work_stream.h>
+
+#include "../tests.h"
+
+
+struct ScratchData
+{};
+
+
+struct CopyData
+{
+ unsigned int computed;
+};
+
+
+struct X
+{
+ void
+ worker(const std::vector<unsigned int>::iterator &i,
+ ScratchData &,
+ CopyData &ad)
+ {
+ ad.computed = *i * 2;
+ }
+
+ void
+ copier(const CopyData &ad)
+ {
+ deallog << ad.computed << std::endl;
+ }
+};
+
+
+void
+test()
+{
+ std::vector<unsigned int> v;
+ for (unsigned int i = 0; i < 20; ++i)
+ v.push_back(i);
+
+ X x;
+ WorkStream::run(v, x, &X::worker, &X::copier, ScratchData(), CopyData());
+}
+
+
+
+int
+main()
+{
+ initlog();
+
+ test();
+}
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2008 - 2017 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.
+//
+// ---------------------------------------------------------------------
+
+
+// test functions in namespace WorkStream
+// This test is based off base/work_stream_01.cc, but for iterator ranges
+
+#include <deal.II/base/work_stream.h>
+
+#include <deal.II/grid/filtered_iterator.h>
+
+#include "../tests.h"
+
+
+struct ScratchData
+{};
+
+
+struct CopyData
+{
+ unsigned int computed;
+};
+
+
+using IteratorType = typename std::vector<unsigned int>::iterator;
+using IteratorRangeType = IteratorRange<IteratorType>;
+
+static_assert(
+ std::is_same<typename IteratorRangeType::iterator, IteratorType>::value,
+ "Iterator types not the same");
+
+struct X
+{
+ void
+ worker(const IteratorType &i, ScratchData &, CopyData &ad)
+ {
+ ad.computed = *i * 2;
+ }
+
+ void
+ copier(const CopyData &ad)
+ {
+ deallog << ad.computed << std::endl;
+ }
+};
+
+
+void
+test()
+{
+ std::vector<unsigned int> v;
+ for (unsigned int i = 0; i < 20; ++i)
+ v.push_back(i);
+
+ const IteratorRangeType iterator_range(v.begin(), v.end());
+
+ X x;
+ WorkStream::run(
+ iterator_range, x, &X::worker, &X::copier, ScratchData(), CopyData());
+}
+
+
+
+int
+main()
+{
+ initlog();
+
+ test();
+}