Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
work_stream.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2019 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_work_stream_h
17 # define dealii_work_stream_h
18 
19 
20 # include <deal.II/base/config.h>
21 
22 # include <deal.II/base/graph_coloring.h>
23 # include <deal.II/base/iterator_range.h>
24 # include <deal.II/base/multithread_info.h>
25 # include <deal.II/base/parallel.h>
26 # include <deal.II/base/template_constraints.h>
27 # include <deal.II/base/thread_local_storage.h>
28 # include <deal.II/base/thread_management.h>
29 
30 # ifdef DEAL_II_WITH_THREADS
31 # include <deal.II/base/thread_management.h>
32 
33 # include <tbb/pipeline.h>
34 # endif
35 
36 # include <functional>
37 # include <iterator>
38 # include <memory>
39 # include <utility>
40 # include <vector>
41 
42 DEAL_II_NAMESPACE_OPEN
43 
44 
45 
149 namespace WorkStream
150 {
151 # ifdef DEAL_II_WITH_THREADS
152 
153  namespace internal
154  {
155  // TODO: The following classes all use std::shared_ptr, but the
156  // correct pointer class would actually be std::unique_ptr. make this
157  // replacement whenever we have a class that provides these semantics
158  // and that is available also as a fall-back whenever via boost or similar
159 
172  namespace Implementation2
173  {
177  template <typename Iterator, typename ScratchData, typename CopyData>
178  class IteratorRangeToItemStream : public tbb::filter
179  {
180  public:
187  struct ItemType
188  {
195  {
196  std::shared_ptr<ScratchData> scratch_data;
197  bool currently_in_use;
198 
203  : currently_in_use(false)
204  {}
205 
206  ScratchDataObject(ScratchData *p, const bool in_use)
207  : scratch_data(p)
208  , currently_in_use(in_use)
209  {}
210 
211  // TODO: when we push back an object to the list of scratch objects,
212  // in
213  // Worker::operator(), we first create an object and then copy
214  // it to the end of this list. this involves having two objects
215  // of the current type having pointers to it, each with their
216  // own currently_in_use flag. there is probably little harm in
217  // this because the original one goes out of scope right away
218  // again, but it's certainly awkward. one way to avoid this
219  // would be to use unique_ptr but we'd need to figure out a way
220  // to use it in non-C++11 mode
222  : scratch_data(o.scratch_data)
223  , currently_in_use(o.currently_in_use)
224  {}
225  };
226 
227 
232  using ScratchDataList = std::list<ScratchDataObject>;
233 
238  std::vector<Iterator> work_items;
239 
245  std::vector<CopyData> copy_datas;
246 
252  unsigned int n_items;
253 
286 
291  const ScratchData *sample_scratch_data;
292 
298 
299 
305  : n_items(0)
306  , scratch_data(nullptr)
307  , sample_scratch_data(nullptr)
308  , currently_in_use(false)
309  {}
310  };
311 
312 
318  IteratorRangeToItemStream(const Iterator & begin,
319  const Iterator & end,
320  const unsigned int buffer_size,
321  const unsigned int chunk_size,
322  const ScratchData &sample_scratch_data,
323  const CopyData & sample_copy_data)
324  : tbb::filter(/*is_serial=*/true)
325  , remaining_iterator_range(begin, end)
326  , item_buffer(buffer_size)
329  {
330  // initialize the elements of the ring buffer
331  for (unsigned int element = 0; element < item_buffer.size();
332  ++element)
333  {
334  Assert(item_buffer[element].n_items == 0, ExcInternalError());
335 
336  item_buffer[element].work_items.resize(
338  item_buffer[element].scratch_data = &thread_local_scratch;
339  item_buffer[element].sample_scratch_data = &sample_scratch_data;
340  item_buffer[element].copy_datas.resize(chunk_size,
341  sample_copy_data);
342  item_buffer[element].currently_in_use = false;
343  }
344  }
345 
346 
350  virtual void *
351  operator()(void *) override
352  {
353  // find first unused item. we know that there must be one
354  // because we have set the maximal number of tokens in flight
355  // and have set the ring buffer to have exactly this size. so
356  // if this function is called, we know that less than the
357  // maximal number of items in currently in flight
358  //
359  // note that we need not lock access to this array since
360  // the current stage is run sequentially and we can therefore
361  // enter the following block only once at any given time.
362  // thus, there can be no race condition between checking that
363  // a flag is false and setting it to true. (there may be
364  // another thread where we release items and set 'false'
365  // flags to 'true', but that too does not produce any
366  // problems)
367  ItemType *current_item = nullptr;
368  for (unsigned int i = 0; i < item_buffer.size(); ++i)
369  if (item_buffer[i].currently_in_use == false)
370  {
371  item_buffer[i].currently_in_use = true;
372  current_item = &item_buffer[i];
373  break;
374  }
375  Assert(current_item != nullptr,
376  ExcMessage("This can't be. There must be a free item!"));
377 
378  // initialize the next item. it may
379  // consist of at most chunk_size
380  // elements
381  current_item->n_items = 0;
382  while ((remaining_iterator_range.first !=
383  remaining_iterator_range.second) &&
384  (current_item->n_items < chunk_size))
385  {
386  current_item->work_items[current_item->n_items] =
388 
389  ++remaining_iterator_range.first;
390  ++current_item->n_items;
391  }
392 
393  if (current_item->n_items == 0)
394  // there were no items
395  // left. terminate the pipeline
396  return nullptr;
397  else
398  return current_item;
399  }
400 
401  private:
406  std::pair<Iterator, Iterator> remaining_iterator_range;
407 
411  std::vector<ItemType> item_buffer;
412 
445 
451  const ScratchData &sample_scratch_data;
452 
459  const unsigned int chunk_size;
460  };
461 
462 
463 
469  template <typename Iterator, typename ScratchData, typename CopyData>
470  class Worker : public tbb::filter
471  {
472  public:
479  const std::function<void(const Iterator &, ScratchData &, CopyData &)>
480  & worker,
481  bool copier_exist = true)
482  : tbb::filter(/* is_serial= */ false)
483  , worker(worker)
485  {}
486 
487 
491  void *
492  operator()(void *item) override
493  {
494  // first unpack the current item
495  using ItemType =
496  typename IteratorRangeToItemStream<Iterator,
497  ScratchData,
498  CopyData>::ItemType;
499 
500  ItemType *current_item = static_cast<ItemType *>(item);
501 
502  // we need to find an unused scratch data object in the list that
503  // corresponds to the current thread and then mark it as used. if
504  // we can't find one, create one
505  //
506  // as discussed in the discussion of the documentation of the
507  // IteratorRangeToItemStream::scratch_data variable, there is no
508  // need to synchronize access to this variable using a mutex
509  // as long as we have no yield-point in between. this means that
510  // we can't take an iterator into the list now and expect it to
511  // still be valid after calling the worker, but we at least do
512  // not have to lock the following section
513  ScratchData *scratch_data = nullptr;
514  {
515  typename ItemType::ScratchDataList &scratch_data_list =
516  current_item->scratch_data->get();
517 
518  // see if there is an unused object. if so, grab it and mark
519  // it as used
520  for (typename ItemType::ScratchDataList::iterator p =
521  scratch_data_list.begin();
522  p != scratch_data_list.end();
523  ++p)
524  if (p->currently_in_use == false)
525  {
526  scratch_data = p->scratch_data.get();
527  p->currently_in_use = true;
528  break;
529  }
530 
531  // if no object was found, create one and mark it as used
532  if (scratch_data == nullptr)
533  {
534  scratch_data =
535  new ScratchData(*current_item->sample_scratch_data);
536 
537  typename ItemType::ScratchDataList::value_type
538  new_scratch_object(scratch_data, true);
539  scratch_data_list.push_back(new_scratch_object);
540  }
541  }
542 
543  // then call the worker function on each element of the chunk we were
544  // given. since these worker functions are called on separate threads,
545  // nothing good can happen if they throw an exception and we are best
546  // off catching it and showing an error message
547  for (unsigned int i = 0; i < current_item->n_items; ++i)
548  {
549  try
550  {
551  if (worker)
552  worker(current_item->work_items[i],
553  *scratch_data,
554  current_item->copy_datas[i]);
555  }
556  catch (const std::exception &exc)
557  {
558  Threads::internal::handle_std_exception(exc);
559  }
560  catch (...)
561  {
562  Threads::internal::handle_unknown_exception();
563  }
564  }
565 
566  // finally mark the scratch object as unused again. as above, there
567  // is no need to lock anything here since the object we work on
568  // is thread-local
569  {
570  typename ItemType::ScratchDataList &scratch_data_list =
571  current_item->scratch_data->get();
572 
573  for (typename ItemType::ScratchDataList::iterator p =
574  scratch_data_list.begin();
575  p != scratch_data_list.end();
576  ++p)
577  if (p->scratch_data.get() == scratch_data)
578  {
579  Assert(p->currently_in_use == true, ExcInternalError());
580  p->currently_in_use = false;
581  }
582  }
583 
584  // if there is no copier, mark current item as usable again
585  if (copier_exist == false)
586  current_item->currently_in_use = false;
587 
588 
589  // then return the original pointer
590  // to the now modified object
591  return item;
592  }
593 
594 
595  private:
600  const std::function<void(const Iterator &, ScratchData &, CopyData &)>
602 
608  };
609 
610 
611 
617  template <typename Iterator, typename ScratchData, typename CopyData>
618  class Copier : public tbb::filter
619  {
620  public:
627  Copier(const std::function<void(const CopyData &)> &copier)
628  : tbb::filter(/*is_serial=*/true)
629  , copier(copier)
630  {}
631 
632 
636  void *
637  operator()(void *item) override
638  {
639  // first unpack the current item
640  using ItemType =
641  typename IteratorRangeToItemStream<Iterator,
642  ScratchData,
643  CopyData>::ItemType;
644 
645  ItemType *current_item = static_cast<ItemType *>(item);
646 
647  // initiate copying data. for the same reasons as in the worker class
648  // above, catch exceptions rather than letting it propagate into
649  // unknown territories
650  for (unsigned int i = 0; i < current_item->n_items; ++i)
651  {
652  try
653  {
654  if (copier)
655  copier(current_item->copy_datas[i]);
656  }
657  catch (const std::exception &exc)
658  {
659  Threads::internal::handle_std_exception(exc);
660  }
661  catch (...)
662  {
663  Threads::internal::handle_unknown_exception();
664  }
665  }
666 
667  // mark current item as usable again
668  current_item->currently_in_use = false;
669 
670 
671  // return an invalid item since we are at the end of the
672  // pipeline
673  return nullptr;
674  }
675 
676 
677  private:
681  const std::function<void(const CopyData &)> copier;
682  };
683 
684  } // namespace Implementation2
685 
686 
694  namespace Implementation3
695  {
701  template <typename Iterator, typename ScratchData, typename CopyData>
703  {
704  std::shared_ptr<ScratchData> scratch_data;
705  std::shared_ptr<CopyData> copy_data;
706  bool currently_in_use;
707 
712  : currently_in_use(false)
713  {}
714 
715  ScratchAndCopyDataObjects(ScratchData *p,
716  CopyData * q,
717  const bool in_use)
718  : scratch_data(p)
719  , copy_data(q)
720  , currently_in_use(in_use)
721  {}
722 
723  // TODO: when we push back an object to the list of scratch objects, in
724  // Worker::operator(), we first create an object and then copy
725  // it to the end of this list. this involves having two objects
726  // of the current type having pointers to it, each with their own
727  // currently_in_use flag. there is probably little harm in this
728  // because the original one goes out of scope right away again, but
729  // it's certainly awkward. one way to avoid this would be to use
730  // unique_ptr but we'd need to figure out a way to use it in
731  // non-C++11 mode
733  : scratch_data(o.scratch_data)
734  , copy_data(o.copy_data)
735  , currently_in_use(o.currently_in_use)
736  {}
737  };
738 
739 
740 
746  template <typename Iterator, typename ScratchData, typename CopyData>
748  {
749  public:
754  const std::function<void(const Iterator &, ScratchData &, CopyData &)>
755  & worker,
756  const std::function<void(const CopyData &)> &copier,
757  const ScratchData & sample_scratch_data,
758  const CopyData & sample_copy_data)
759  : worker(worker)
760  , copier(copier)
762  , sample_copy_data(sample_copy_data)
763  {}
764 
765 
770  void
771  operator()(const tbb::blocked_range<
772  typename std::vector<Iterator>::const_iterator> &range)
773  {
774  // we need to find an unused scratch and corresponding copy
775  // data object in the list that corresponds to the current
776  // thread and then mark it as used. If we can't find one,
777  // create one as discussed in the discussion of the documentation
778  // of the IteratorRangeToItemStream::scratch_data variable,
779  // there is no need to synchronize access to this variable
780  // using a mutex as long as we have no yield-point in between.
781  // This means that we can't take an iterator into the list
782  // now and expect it to still be valid after calling the worker,
783  // but we at least do not have to lock the following section.
784  ScratchData *scratch_data = nullptr;
785  CopyData * copy_data = nullptr;
786  {
787  ScratchAndCopyDataList &scratch_and_copy_data_list = data.get();
788 
789  // see if there is an unused object. if so, grab it and mark
790  // it as used
791  for (typename ScratchAndCopyDataList::iterator p =
792  scratch_and_copy_data_list.begin();
793  p != scratch_and_copy_data_list.end();
794  ++p)
795  if (p->currently_in_use == false)
796  {
797  scratch_data = p->scratch_data.get();
798  copy_data = p->copy_data.get();
799  p->currently_in_use = true;
800  break;
801  }
802 
803  // if no element in the list was found, create one and mark it as
804  // used
805  if (scratch_data == nullptr)
806  {
807  Assert(copy_data == nullptr, ExcInternalError());
808  scratch_data = new ScratchData(sample_scratch_data);
809  copy_data = new CopyData(sample_copy_data);
810 
811  scratch_and_copy_data_list.emplace_back(scratch_data,
812  copy_data,
813  true);
814  }
815  }
816 
817  // then call the worker and copier functions on each
818  // element of the chunk we were given.
819  for (typename std::vector<Iterator>::const_iterator p = range.begin();
820  p != range.end();
821  ++p)
822  {
823  try
824  {
825  if (worker)
826  worker(*p, *scratch_data, *copy_data);
827  if (copier)
828  copier(*copy_data);
829  }
830  catch (const std::exception &exc)
831  {
832  Threads::internal::handle_std_exception(exc);
833  }
834  catch (...)
835  {
836  Threads::internal::handle_unknown_exception();
837  }
838  }
839 
840  // finally mark the scratch object as unused again. as above, there
841  // is no need to lock anything here since the object we work on
842  // is thread-local
843  {
844  ScratchAndCopyDataList &scratch_and_copy_data_list = data.get();
845 
846  for (typename ScratchAndCopyDataList::iterator p =
847  scratch_and_copy_data_list.begin();
848  p != scratch_and_copy_data_list.end();
849  ++p)
850  if (p->scratch_data.get() == scratch_data)
851  {
852  Assert(p->currently_in_use == true, ExcInternalError());
853  p->currently_in_use = false;
854  }
855  }
856  }
857 
858  private:
859  using ScratchAndCopyDataObjects = typename Implementation3::
860  ScratchAndCopyDataObjects<Iterator, ScratchData, CopyData>;
861 
866  using ScratchAndCopyDataList = std::list<ScratchAndCopyDataObjects>;
867 
869 
874  const std::function<void(const Iterator &, ScratchData &, CopyData &)>
876 
881  const std::function<void(const CopyData &)> copier;
882 
886  const ScratchData &sample_scratch_data;
887  const CopyData & sample_copy_data;
888  };
889  } // namespace Implementation3
890 
891  } // namespace internal
892 
893 
894 # endif // DEAL_II_WITH_THREADS
895 
896 
931  template <typename Worker,
932  typename Copier,
933  typename Iterator,
934  typename ScratchData,
935  typename CopyData>
936  void
937  run(const std::vector<std::vector<Iterator>> &colored_iterators,
938  Worker worker,
939  Copier copier,
940  const ScratchData & sample_scratch_data,
941  const CopyData & sample_copy_data,
942  const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
943  const unsigned int chunk_size = 8);
944 
945 
980  template <typename Worker,
981  typename Copier,
982  typename Iterator,
983  typename ScratchData,
984  typename CopyData>
985  void
986  run(const Iterator & begin,
987  const typename identity<Iterator>::type &end,
988  Worker worker,
989  Copier copier,
990  const ScratchData & sample_scratch_data,
991  const CopyData & sample_copy_data,
992  const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
993  const unsigned int chunk_size = 8)
994  {
995  Assert(queue_length > 0,
996  ExcMessage("The queue length must be at least one, and preferably "
997  "larger than the number of processors on this system."));
998  (void)queue_length; // removes -Wunused-parameter warning in optimized mode
999  Assert(chunk_size > 0, ExcMessage("The chunk_size must be at least one."));
1000  (void)chunk_size; // removes -Wunused-parameter warning in optimized mode
1001 
1002  // if no work then skip. (only use operator!= for iterators since we may
1003  // not have an equality comparison operator)
1004  if (!(begin != end))
1005  return;
1006 
1007  // we want to use TBB if we have support and if it is not disabled at
1008  // runtime:
1009 # ifdef DEAL_II_WITH_THREADS
1010  if (MultithreadInfo::n_threads() == 1)
1011 # endif
1012  {
1013  // need to copy the sample since it is marked const
1014  ScratchData scratch_data = sample_scratch_data;
1015  CopyData copy_data = sample_copy_data; // NOLINT
1016 
1017  for (Iterator i = begin; i != end; ++i)
1018  {
1019  // need to check if the function is not the zero function. To
1020  // check zero-ness, create a C++ function out of it and check that
1021  if (static_cast<const std::function<
1022  void(const Iterator &, ScratchData &, CopyData &)> &>(worker))
1023  worker(i, scratch_data, copy_data);
1024  if (static_cast<const std::function<void(const CopyData &)> &>(
1025  copier))
1026  copier(copy_data);
1027  }
1028  }
1029 # ifdef DEAL_II_WITH_THREADS
1030  else // have TBB and use more than one thread
1031  {
1032  // Check that the copier exist
1033  if (static_cast<const std::function<void(const CopyData &)> &>(copier))
1034  {
1035  // create the three stages of the pipeline
1036  internal::Implementation2::
1037  IteratorRangeToItemStream<Iterator, ScratchData, CopyData>
1038  iterator_range_to_item_stream(begin,
1039  end,
1040  queue_length,
1041  chunk_size,
1042  sample_scratch_data,
1043  sample_copy_data);
1044 
1046  worker_filter(worker);
1048  copier_filter(copier);
1049 
1050  // now create a pipeline from these stages
1051  tbb::pipeline assembly_line;
1052  assembly_line.add_filter(iterator_range_to_item_stream);
1053  assembly_line.add_filter(worker_filter);
1054  assembly_line.add_filter(copier_filter);
1055 
1056  // and run it
1057  assembly_line.run(queue_length);
1058 
1059  assembly_line.clear();
1060  }
1061  else
1062  {
1063  // there is no copier function. in this case, we have an
1064  // embarrassingly parallel problem where we can
1065  // essentially apply parallel_for. because parallel_for
1066  // requires subdividing the range for which operator- is
1067  // necessary between iterators, it is often inefficient to
1068  // apply it directly to cell ranges and similar iterator
1069  // types for which operator- is expensive or, in fact,
1070  // nonexistent. rather, in that case, we simply copy the
1071  // iterators into a large array and use operator- on
1072  // iterators to this array of iterators.
1073  //
1074  // instead of duplicating code, this is essentially the
1075  // same situation we have in Implementation3 below, so we
1076  // just defer to that place
1077  std::vector<std::vector<Iterator>> all_iterators(1);
1078  for (Iterator p = begin; p != end; ++p)
1079  all_iterators[0].push_back(p);
1080 
1081  run(all_iterators,
1082  worker,
1083  copier,
1084  sample_scratch_data,
1085  sample_copy_data,
1086  queue_length,
1087  chunk_size);
1088  }
1089  }
1090 # endif
1091  }
1092 
1093 
1094 
1102  template <typename Worker,
1103  typename Copier,
1104  typename IteratorRangeType,
1105  typename ScratchData,
1106  typename CopyData,
1107  typename = typename std::enable_if<
1109  void
1110  run(IteratorRangeType iterator_range,
1111  Worker worker,
1112  Copier copier,
1113  const ScratchData &sample_scratch_data,
1114  const CopyData & sample_copy_data,
1115  const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
1116  const unsigned int chunk_size = 8)
1117  {
1118  // Call the function above
1119  run(iterator_range.begin(),
1120  iterator_range.end(),
1121  worker,
1122  copier,
1123  sample_scratch_data,
1124  sample_copy_data,
1125  queue_length,
1126  chunk_size);
1127  }
1128 
1129 
1130 
1134  template <typename Worker,
1135  typename Copier,
1136  typename Iterator,
1137  typename ScratchData,
1138  typename CopyData>
1139  void
1140  run(const IteratorRange<Iterator> &iterator_range,
1141  Worker worker,
1142  Copier copier,
1143  const ScratchData & sample_scratch_data,
1144  const CopyData & sample_copy_data,
1145  const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
1146  const unsigned int chunk_size = 8)
1147  {
1148  // Call the function above
1149  run(iterator_range.begin(),
1150  iterator_range.end(),
1151  worker,
1152  copier,
1153  sample_scratch_data,
1154  sample_copy_data,
1155  queue_length,
1156  chunk_size);
1157  }
1158 
1159 
1160  // Implementation 3:
1161  template <typename Worker,
1162  typename Copier,
1163  typename Iterator,
1164  typename ScratchData,
1165  typename CopyData>
1166  void
1167  run(const std::vector<std::vector<Iterator>> &colored_iterators,
1168  Worker worker,
1169  Copier copier,
1170  const ScratchData & sample_scratch_data,
1171  const CopyData & sample_copy_data,
1172  const unsigned int queue_length,
1173  const unsigned int chunk_size)
1174  {
1175  Assert(queue_length > 0,
1176  ExcMessage("The queue length must be at least one, and preferably "
1177  "larger than the number of processors on this system."));
1178  (void)queue_length; // removes -Wunused-parameter warning in optimized mode
1179  Assert(chunk_size > 0, ExcMessage("The chunk_size must be at least one."));
1180  (void)chunk_size; // removes -Wunused-parameter warning in optimized mode
1181 
1182  // we want to use TBB if we have support and if it is not disabled at
1183  // runtime:
1184 # ifdef DEAL_II_WITH_THREADS
1185  if (MultithreadInfo::n_threads() == 1)
1186 # endif
1187  {
1188  // need to copy the sample since it is marked const
1189  ScratchData scratch_data = sample_scratch_data;
1190  CopyData copy_data = sample_copy_data; // NOLINT
1191 
1192  for (unsigned int color = 0; color < colored_iterators.size(); ++color)
1193  for (typename std::vector<Iterator>::const_iterator p =
1194  colored_iterators[color].begin();
1195  p != colored_iterators[color].end();
1196  ++p)
1197  {
1198  // need to check if the function is not the zero function. To
1199  // check zero-ness, create a C++ function out of it and check that
1200  if (static_cast<const std::function<void(
1201  const Iterator &, ScratchData &, CopyData &)> &>(worker))
1202  worker(*p, scratch_data, copy_data);
1203  if (static_cast<const std::function<void(const CopyData &)> &>(
1204  copier))
1205  copier(copy_data);
1206  }
1207  }
1208 # ifdef DEAL_II_WITH_THREADS
1209  else // have TBB and use more than one thread
1210  {
1211  // loop over the various colors of what we're given
1212  for (unsigned int color = 0; color < colored_iterators.size(); ++color)
1213  if (colored_iterators[color].size() > 0)
1214  {
1215  using WorkerAndCopier = internal::Implementation3::
1216  WorkerAndCopier<Iterator, ScratchData, CopyData>;
1217 
1218  using RangeType = typename std::vector<Iterator>::const_iterator;
1219 
1220  WorkerAndCopier worker_and_copier(worker,
1221  copier,
1222  sample_scratch_data,
1223  sample_copy_data);
1224 
1226  colored_iterators[color].begin(),
1227  colored_iterators[color].end(),
1228  std::bind(&WorkerAndCopier::operator(),
1229  std::ref(worker_and_copier),
1230  std::placeholders::_1),
1231  chunk_size);
1232  }
1233  }
1234 # endif
1235  }
1236 
1237 
1238 
1268  template <typename MainClass,
1269  typename Iterator,
1270  typename ScratchData,
1271  typename CopyData>
1272  void
1273  run(const Iterator & begin,
1274  const typename identity<Iterator>::type &end,
1275  MainClass & main_object,
1276  void (MainClass::*worker)(const Iterator &, ScratchData &, CopyData &),
1277  void (MainClass::*copier)(const CopyData &),
1278  const ScratchData &sample_scratch_data,
1279  const CopyData & sample_copy_data,
1280  const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
1281  const unsigned int chunk_size = 8)
1282  {
1283  // forward to the other function
1284  run(begin,
1285  end,
1286  std::bind(worker,
1287  std::ref(main_object),
1288  std::placeholders::_1,
1289  std::placeholders::_2,
1290  std::placeholders::_3),
1291  std::bind(copier, std::ref(main_object), std::placeholders::_1),
1292  sample_scratch_data,
1293  sample_copy_data,
1294  queue_length,
1295  chunk_size);
1296  }
1297 
1298 
1299  template <typename MainClass,
1300  typename Iterator,
1301  typename ScratchData,
1302  typename CopyData>
1303  void
1304  run(const IteratorOverIterators<Iterator> & begin,
1305  const IteratorOverIterators<typename identity<Iterator>::type> &end,
1306  MainClass &main_object,
1307  void (MainClass::*worker)(const Iterator &, ScratchData &, CopyData &),
1308  void (MainClass::*copier)(const CopyData &),
1309  const ScratchData &sample_scratch_data,
1310  const CopyData & sample_copy_data,
1311  const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
1312  const unsigned int chunk_size = 8)
1313  {
1314  // forward to the other function
1315  run(begin,
1316  end,
1317  std::bind(worker,
1318  std::ref(main_object),
1319  std::placeholders::_1,
1320  std::placeholders::_2,
1321  std::placeholders::_3),
1322  std::bind(copier, std::ref(main_object), std::placeholders::_1),
1323  sample_scratch_data,
1324  sample_copy_data,
1325  queue_length,
1326  chunk_size);
1327  }
1328 
1329 
1330 
1338  template <typename MainClass,
1339  typename IteratorRangeType,
1340  typename ScratchData,
1341  typename CopyData,
1342  typename = typename std::enable_if<
1344  void
1345  run(IteratorRangeType iterator_range,
1346  MainClass & main_object,
1347  void (MainClass::*worker)(
1349  ScratchData &,
1350  CopyData &),
1351  void (MainClass::*copier)(const CopyData &),
1352  const ScratchData &sample_scratch_data,
1353  const CopyData & sample_copy_data,
1354  const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
1355  const unsigned int chunk_size = 8)
1356  {
1357  // Call the function above
1358  run(std::begin(iterator_range),
1359  std::end(iterator_range),
1360  main_object,
1361  worker,
1362  copier,
1363  sample_scratch_data,
1364  sample_copy_data,
1365  queue_length,
1366  chunk_size);
1367  }
1368 
1369 
1370 
1374  template <typename MainClass,
1375  typename Iterator,
1376  typename ScratchData,
1377  typename CopyData>
1378  void
1380  MainClass & main_object,
1381  void (MainClass::*worker)(const Iterator &, ScratchData &, CopyData &),
1382  void (MainClass::*copier)(const CopyData &),
1383  const ScratchData &sample_scratch_data,
1384  const CopyData & sample_copy_data,
1385  const unsigned int queue_length = 2 * MultithreadInfo::n_threads(),
1386  const unsigned int chunk_size = 8)
1387  {
1388  // Call the function above
1389  run(std::begin(iterator_range),
1390  std::end(iterator_range),
1391  main_object,
1392  worker,
1393  copier,
1394  sample_scratch_data,
1395  sample_copy_data,
1396  queue_length,
1397  chunk_size);
1398  }
1399 
1400 } // namespace WorkStream
1401 
1402 
1403 
1404 DEAL_II_NAMESPACE_CLOSE
1405 
1406 
1407 
1408 //---------------------------- work_stream.h ---------------------------
1409 // end of #ifndef dealii_work_stream_h
1410 #endif
1411 //---------------------------- work_stream.h ---------------------------
void * operator()(void *item) override
Definition: work_stream.h:637
const std::function< void(const Iterator &, ScratchData &, CopyData &)> worker
Definition: work_stream.h:601
const std::function< void(const CopyData &)> copier
Definition: work_stream.h:881
std::list< ScratchAndCopyDataObjects > ScratchAndCopyDataList
Definition: work_stream.h:866
IteratorOverIterators begin()
void parallel_for(Iterator x_begin, Iterator x_end, const Functor &functor, const unsigned int grainsize)
Definition: parallel.h:156
const std::function< void(const Iterator &, ScratchData &, CopyData &)> worker
Definition: work_stream.h:875
void * operator()(void *item) override
Definition: work_stream.h:492
void operator()(const tbb::blocked_range< typename std::vector< Iterator >::const_iterator > &range)
Definition: work_stream.h:771
Threads::ThreadLocalStorage< ScratchDataList > * scratch_data
Definition: work_stream.h:285
IteratorRangeToItemStream(const Iterator &begin, const Iterator &end, const unsigned int buffer_size, const unsigned int chunk_size, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data)
Definition: work_stream.h:318
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1407
Worker(const std::function< void(const Iterator &, ScratchData &, CopyData &)> &worker, bool copier_exist=true)
Definition: work_stream.h:478
WorkerAndCopier(const std::function< void(const Iterator &, ScratchData &, CopyData &)> &worker, const std::function< void(const CopyData &)> &copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data)
Definition: work_stream.h:753
const std::function< void(const CopyData &)> copier
Definition: work_stream.h:681
void run(const std::vector< std::vector< Iterator >> &colored_iterators, 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)
Definition: work_stream.h:1167
static unsigned int n_threads()
Threads::ThreadLocalStorage< typename ItemType::ScratchDataList > thread_local_scratch
Definition: work_stream.h:444
Copier(const std::function< void(const CopyData &)> &copier)
Definition: work_stream.h:627
static ::ExceptionBase & ExcInternalError()
IteratorOverIterators end() const