Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
parallel.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_parallel_h
17 #define dealii_parallel_h
18 
19 
20 #include <deal.II/base/config.h>
21 
26 
27 #include <cstddef>
28 #include <functional>
29 #include <memory>
30 #include <tuple>
31 
32 #ifdef DEAL_II_WITH_THREADS
33 # include <tbb/blocked_range.h>
34 # include <tbb/parallel_for.h>
35 # include <tbb/parallel_reduce.h>
36 # include <tbb/partitioner.h>
37 #endif
38 
39 
40 // TODO[WB]: allow calling functions to pass along a tbb::affinity_partitioner
41 // object to ensure that subsequent calls use the same cache lines
42 
44 
45 namespace parallel
46 {
47  namespace internal
48  {
53  template <typename Number>
55  {
56  static const bool value = true;
57  };
58 
59 #ifdef __INTEL_COMPILER
60  // Disable long double SIMD instructions on ICC. This is to work around a
61  // bug that generates wrong code at least up to intel 15 (see
62  // tests/lac/vector-vector, tests/lac/intel-15-bug, and the discussion at
63  // https://github.com/dealii/dealii/issues/598).
64  template <>
65  struct EnableOpenMPSimdFor<long double>
66  {
67  static const bool value = false;
68  };
69 #endif
70 
71 
72 
77  template <typename F>
78  struct Body
79  {
83  Body(const F &f)
84  : f(f)
85  {}
86 
87  template <typename Range>
88  void
89  operator()(const Range &range) const
90  {
91  for (typename Range::const_iterator p = range.begin(); p != range.end();
92  ++p)
93  apply(f, *p);
94  }
95 
96  private:
100  const F f;
101 
105  template <typename I1, typename I2>
106  static void
107  apply(const F &f, const std::tuple<I1, I2> &p)
108  {
109  *std::get<1>(p) = f(*std::get<0>(p));
110  }
111 
115  template <typename I1, typename I2, typename I3>
116  static void
117  apply(const F &f, const std::tuple<I1, I2, I3> &p)
118  {
119  *std::get<2>(p) = f(*std::get<0>(p), *std::get<1>(p));
120  }
121 
125  template <typename I1, typename I2, typename I3, typename I4>
126  static void
127  apply(const F &f, const std::tuple<I1, I2, I3, I4> &p)
128  {
129  *std::get<3>(p) = f(*std::get<0>(p), *std::get<1>(p), *std::get<2>(p));
130  }
131  };
132 
133 
134 
141  template <typename F>
142  Body<F>
143  make_body(const F &f)
144  {
145  return Body<F>(f);
146  }
147 
148 
149 
150 #ifdef DEAL_II_WITH_THREADS
151 
154  template <typename Iterator, typename Functor>
155  void
157  Iterator x_end,
158  const Functor & functor,
159  const unsigned int grainsize)
160  {
161  tbb::parallel_for(tbb::blocked_range<Iterator>(x_begin, x_end, grainsize),
162  functor,
163  tbb::auto_partitioner());
164  }
165 
166 
167 
171  template <typename Iterator, typename Functor>
172  void
174  Iterator x_end,
175  const Functor & functor,
176  const unsigned int grainsize,
177  const std::shared_ptr<tbb::affinity_partitioner> &partitioner)
178  {
179  tbb::parallel_for(tbb::blocked_range<Iterator>(x_begin, x_end, grainsize),
180  functor,
181  *partitioner);
182  }
183 #endif
184  } // namespace internal
185 
209  template <typename InputIterator, typename OutputIterator, typename Predicate>
210  void
211  transform(const InputIterator &begin_in,
212  const InputIterator &end_in,
213  OutputIterator out,
214  const Predicate & predicate,
215  const unsigned int grainsize)
216  {
217 #ifndef DEAL_II_WITH_THREADS
218  // make sure we don't get compiler
219  // warnings about unused arguments
220  (void)grainsize;
221 
222  for (OutputIterator in = begin_in; in != end_in;)
223  *out++ = predicate(*in++);
224 #else
225  using Iterators = std::tuple<InputIterator, OutputIterator>;
226  using SyncIterators = SynchronousIterators<Iterators>;
227  Iterators x_begin(begin_in, out);
228  Iterators x_end(end_in, OutputIterator());
229  internal::parallel_for(SyncIterators(x_begin),
230  SyncIterators(x_end),
231  internal::make_body(predicate),
232  grainsize);
233 #endif
234  }
235 
236 
237 
261  template <typename InputIterator1,
262  typename InputIterator2,
263  typename OutputIterator,
264  typename Predicate>
265  void
266  transform(const InputIterator1 &begin_in1,
267  const InputIterator1 &end_in1,
268  InputIterator2 in2,
269  OutputIterator out,
270  const Predicate & predicate,
271  const unsigned int grainsize)
272  {
273 #ifndef DEAL_II_WITH_THREADS
274  // make sure we don't get compiler
275  // warnings about unused arguments
276  (void)grainsize;
277 
278  for (OutputIterator in1 = begin_in1; in1 != end_in1;)
279  *out++ = predicate(*in1++, *in2++);
280 #else
281  using Iterators =
282  std::tuple<InputIterator1, InputIterator2, OutputIterator>;
283  using SyncIterators = SynchronousIterators<Iterators>;
284  Iterators x_begin(begin_in1, in2, out);
285  Iterators x_end(end_in1, InputIterator2(), OutputIterator());
286  internal::parallel_for(SyncIterators(x_begin),
287  SyncIterators(x_end),
288  internal::make_body(predicate),
289  grainsize);
290 #endif
291  }
292 
293 
294 
318  template <typename InputIterator1,
319  typename InputIterator2,
320  typename InputIterator3,
321  typename OutputIterator,
322  typename Predicate>
323  void
324  transform(const InputIterator1 &begin_in1,
325  const InputIterator1 &end_in1,
326  InputIterator2 in2,
327  InputIterator3 in3,
328  OutputIterator out,
329  const Predicate & predicate,
330  const unsigned int grainsize)
331  {
332 #ifndef DEAL_II_WITH_THREADS
333  // make sure we don't get compiler
334  // warnings about unused arguments
335  (void)grainsize;
336 
337  for (OutputIterator in1 = begin_in1; in1 != end_in1;)
338  *out++ = predicate(*in1++, *in2++, *in3++);
339 #else
340  using Iterators = std::
341  tuple<InputIterator1, InputIterator2, InputIterator3, OutputIterator>;
342  using SyncIterators = SynchronousIterators<Iterators>;
343  Iterators x_begin(begin_in1, in2, in3, out);
344  Iterators x_end(end_in1,
345  InputIterator2(),
346  InputIterator3(),
347  OutputIterator());
348  internal::parallel_for(SyncIterators(x_begin),
349  SyncIterators(x_end),
350  internal::make_body(predicate),
351  grainsize);
352 #endif
353  }
354 
355 
356  namespace internal
357  {
358 #ifdef DEAL_II_WITH_THREADS
359 
363  template <typename RangeType, typename Function>
364  void
365  apply_to_subranges(const tbb::blocked_range<RangeType> &range,
366  const Function & f)
367  {
368  f(range.begin(), range.end());
369  }
370 #endif
371  } // namespace internal
372 
373 
443  template <typename RangeType, typename Function>
444  void
445  apply_to_subranges(const RangeType & begin,
446  const typename identity<RangeType>::type &end,
447  const Function & f,
448  const unsigned int grainsize)
449  {
450 #ifndef DEAL_II_WITH_THREADS
451  // make sure we don't get compiler
452  // warnings about unused arguments
453  (void)grainsize;
454 
455 # ifndef DEAL_II_BIND_NO_CONST_OP_PARENTHESES
456  f(begin, end);
457 # else
458  // work around a problem with MS VC++ where there is no const
459  // operator() in 'Function' if 'Function' is the result of std::bind
460  Function ff = f;
461  ff(begin, end);
462 # endif
463 #else
465  end,
466  [&f](const tbb::blocked_range<RangeType> &range) {
467  internal::apply_to_subranges<RangeType, Function>(
468  range, f);
469  },
470  grainsize);
471 #endif
472  }
473 
474 
475 
504  {
509  virtual ~ParallelForInteger() = default;
510 
519  void
520  apply_parallel(const std::size_t begin,
521  const std::size_t end,
522  const std::size_t minimum_parallel_grain_size) const;
523 
530  virtual void
531  apply_to_subrange(const std::size_t, const std::size_t) const = 0;
532  };
533 
534 
535 
536  namespace internal
537  {
538 #ifdef DEAL_II_WITH_THREADS
539 
545  template <typename ResultType, typename Function>
547  {
551  ResultType result;
552 
562  template <typename Reductor>
564  const Reductor & reductor,
565  const ResultType neutral_element = ResultType())
567  , f(f)
569  , reductor(reductor)
570  {}
571 
577  , f(r.f)
579  , reductor(r.reductor)
580  {}
581 
586  void
588  {
590  }
591 
595  template <typename RangeType>
596  void
597  operator()(const tbb::blocked_range<RangeType> &range)
598  {
599  result = reductor(result, f(range.begin(), range.end()));
600  }
601 
602  private:
606  const Function f;
607 
613  const ResultType neutral_element;
614 
619  const std::function<ResultType(ResultType, ResultType)> reductor;
620  };
621 #endif
622  } // namespace internal
623 
624 
685  template <typename ResultType, typename RangeType, typename Function>
686  ResultType
688  const RangeType & begin,
689  const typename identity<RangeType>::type &end,
690  const unsigned int grainsize)
691  {
692 #ifndef DEAL_II_WITH_THREADS
693  // make sure we don't get compiler
694  // warnings about unused arguments
695  (void)grainsize;
696 
697 # ifndef DEAL_II_BIND_NO_CONST_OP_PARENTHESES
698  return f(begin, end);
699 # else
700  // work around a problem with MS VC++ where there is no const
701  // operator() in 'Function' if 'Function' is the result of std::bind
702  Function ff = f;
703  return ff(begin, end);
704 # endif
705 #else
707  f, std::plus<ResultType>(), 0);
708  tbb::parallel_reduce(tbb::blocked_range<RangeType>(begin, end, grainsize),
709  reductor,
710  tbb::auto_partitioner());
711  return reductor.result;
712 #endif
713  }
714 
715 
716  // --------------------- for loop affinity partitioner -----------------------
717 
727  namespace internal
728  {
730  {
731  public:
735  TBBPartitioner();
736 
737 #ifdef DEAL_II_WITH_THREADS
738 
742  ~TBBPartitioner();
743 
750  std::shared_ptr<tbb::affinity_partitioner>
752 
758  void
759  release_one_partitioner(std::shared_ptr<tbb::affinity_partitioner> &p);
760 
761  private:
766  std::shared_ptr<tbb::affinity_partitioner> my_partitioner;
767 
772  bool in_use;
773 
777  std::mutex mutex;
778 #endif
779  };
780  } // namespace internal
781 } // namespace parallel
782 
783 
784 namespace internal
785 {
786  namespace VectorImplementation
787  {
802  extern unsigned int minimum_parallel_grain_size;
803  } // namespace VectorImplementation
804 
805 
806  namespace SparseMatrixImplementation
807  {
813  extern unsigned int minimum_parallel_grain_size;
814  } // namespace SparseMatrixImplementation
815 
816 } // end of namespace internal
817 
818 
819 /* --------------------------- inline functions ------------------------- */
820 
821 namespace parallel
822 {
823 #ifdef DEAL_II_WITH_THREADS
824 
825  namespace internal
826  {
832  {
834  : worker_(worker)
835  {}
836 
837  void
838  operator()(const tbb::blocked_range<std::size_t> &range) const
839  {
840  worker_.apply_to_subrange(range.begin(), range.end());
841  }
842 
844  };
845  } // namespace internal
846 
847 #endif
848 
849 
850  inline void
852  const std::size_t begin,
853  const std::size_t end,
854  const std::size_t minimum_parallel_grain_size) const
855  {
856 #ifndef DEAL_II_WITH_THREADS
857  // make sure we don't get compiler
858  // warnings about unused arguments
860 
862 #else
863  internal::ParallelForWrapper worker(*this);
865 #endif
866  }
867 
868 } // end of namespace parallel
869 
871 
872 #endif
identity::type
T type
Definition: template_constraints.h:270
Physics::Elasticity::Kinematics::F
Tensor< 2, dim, Number > F(const Tensor< 2, dim, Number > &Grad_u)
parallel::internal::TBBPartitioner::my_partitioner
std::shared_ptr< tbb::affinity_partitioner > my_partitioner
Definition: parallel.h:766
parallel::internal::TBBPartitioner::in_use
bool in_use
Definition: parallel.h:772
SynchronousIterators
Definition: synchronous_iterator.h:52
parallel::internal::make_body
Body< F > make_body(const F &f)
Definition: parallel.h:143
parallel::internal::Body::apply
static void apply(const F &f, const std::tuple< I1, I2, I3 > &p)
Definition: parallel.h:117
parallel::internal::TBBPartitioner::mutex
std::mutex mutex
Definition: parallel.h:777
parallel::transform
void transform(const InputIterator &begin_in, const InputIterator &end_in, OutputIterator out, const Predicate &predicate, const unsigned int grainsize)
Definition: parallel.h:211
parallel::internal::ParallelForWrapper::ParallelForWrapper
ParallelForWrapper(const parallel::ParallelForInteger &worker)
Definition: parallel.h:833
parallel::internal::Body::Body
Body(const F &f)
Definition: parallel.h:83
parallel::internal::ReductionOnSubranges::ReductionOnSubranges
ReductionOnSubranges(const ReductionOnSubranges &r, tbb::split)
Definition: parallel.h:575
parallel::internal::Body
Definition: parallel.h:78
parallel::ParallelForInteger::apply_to_subrange
virtual void apply_to_subrange(const std::size_t, const std::size_t) const =0
parallel::internal::TBBPartitioner::TBBPartitioner
TBBPartitioner()
Definition: parallel.cc:57
internal::SparseMatrixImplementation::minimum_parallel_grain_size
unsigned int minimum_parallel_grain_size
Definition: parallel.cc:48
parallel::internal::ReductionOnSubranges::operator()
void operator()(const tbb::blocked_range< RangeType > &range)
Definition: parallel.h:597
parallel::internal::ParallelForWrapper::operator()
void operator()(const tbb::blocked_range< std::size_t > &range) const
Definition: parallel.h:838
parallel::internal::ReductionOnSubranges::reductor
const std::function< ResultType(ResultType, ResultType)> reductor
Definition: parallel.h:619
parallel::internal::TBBPartitioner::release_one_partitioner
void release_one_partitioner(std::shared_ptr< tbb::affinity_partitioner > &p)
Definition: parallel.cc:88
thread_management.h
parallel::internal::ReductionOnSubranges::f
const Function f
Definition: parallel.h:606
TransposeTableIterators::Iterator
MatrixTableIterators::Iterator< TransposeTable< T >, Constness, MatrixTableIterators::Storage::column_major > Iterator
Definition: table.h:1913
parallel::internal::Body::apply
static void apply(const F &f, const std::tuple< I1, I2, I3, I4 > &p)
Definition: parallel.h:127
parallel::internal::TBBPartitioner
Definition: parallel.h:729
double
TrilinosWrappers::internal::begin
VectorType::value_type * begin(VectorType &V)
Definition: trilinos_sparse_matrix.cc:51
internal::VectorOperations::parallel_reduce
void parallel_reduce(const Operation &op, const size_type start, const size_type end, ResultType &result, const std::shared_ptr<::parallel::internal::TBBPartitioner > &partitioner)
Definition: vector_operations_internal.h:1353
parallel::internal::Body::f
const F f
Definition: parallel.h:100
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
parallel::internal::ReductionOnSubranges::ReductionOnSubranges
ReductionOnSubranges(const Function &f, const Reductor &reductor, const ResultType neutral_element=ResultType())
Definition: parallel.h:563
parallel::internal::Body::apply
static void apply(const F &f, const std::tuple< I1, I2 > &p)
Definition: parallel.h:107
internal::VectorImplementation::minimum_parallel_grain_size
unsigned int minimum_parallel_grain_size
Definition: parallel.cc:34
parallel::ParallelForInteger::~ParallelForInteger
virtual ~ParallelForInteger()=default
TrilinosWrappers::internal::end
VectorType::value_type * end(VectorType &V)
Definition: trilinos_sparse_matrix.cc:65
parallel::internal::TBBPartitioner::acquire_one_partitioner
std::shared_ptr< tbb::affinity_partitioner > acquire_one_partitioner()
Definition: parallel.cc:75
parallel::internal::Body::operator()
void operator()(const Range &range) const
Definition: parallel.h:89
exceptions.h
parallel::accumulate_from_subranges
ResultType accumulate_from_subranges(const Function &f, const RangeType &begin, const typename identity< RangeType >::type &end, const unsigned int grainsize)
Definition: parallel.h:687
AdaptationStrategies::Refinement::split
std::vector< value_type > split(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)
parallel::internal::ReductionOnSubranges::result
ResultType result
Definition: parallel.h:551
parallel::internal::ParallelForWrapper::worker_
const parallel::ParallelForInteger & worker_
Definition: parallel.h:843
parallel::internal::EnableOpenMPSimdFor
Definition: parallel.h:54
parallel::internal::TBBPartitioner::~TBBPartitioner
~TBBPartitioner()
Definition: parallel.cc:64
parallel::internal::ParallelForWrapper
Definition: parallel.h:831
template_constraints.h
config.h
Function
Definition: function.h:151
internal
Definition: aligned_vector.h:369
parallel::internal::ReductionOnSubranges::neutral_element
const ResultType neutral_element
Definition: parallel.h:613
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
parallel::internal::ReductionOnSubranges
Definition: parallel.h:546
parallel::internal::EnableOpenMPSimdFor::value
static const bool value
Definition: parallel.h:56
parallel
Definition: distributed.h:416
parallel::ParallelForInteger::apply_parallel
void apply_parallel(const std::size_t begin, const std::size_t end, const std::size_t minimum_parallel_grain_size) const
Definition: parallel.h:851
parallel::internal::ReductionOnSubranges::join
void join(const ReductionOnSubranges &r)
Definition: parallel.h:587
parallel::apply_to_subranges
void apply_to_subranges(const RangeType &begin, const typename identity< RangeType >::type &end, const Function &f, const unsigned int grainsize)
Definition: parallel.h:445
parallel::ParallelForInteger
Definition: parallel.h:503
synchronous_iterator.h
parallel::internal::apply_to_subranges
void apply_to_subranges(const tbb::blocked_range< RangeType > &range, const Function &f)
Definition: parallel.h:365
parallel::internal::parallel_for
void parallel_for(Iterator x_begin, Iterator x_end, const Functor &functor, const unsigned int grainsize, const std::shared_ptr< tbb::affinity_partitioner > &partitioner)
Definition: parallel.h:173
parallel::internal::parallel_for
void parallel_for(Iterator x_begin, Iterator x_end, const Functor &functor, const unsigned int grainsize)
Definition: parallel.h:156