Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
task_info.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2018 - 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 
17 #include <deal.II/base/conditional_ostream.h>
18 #include <deal.II/base/memory_consumption.h>
19 #include <deal.II/base/mpi.h>
20 #include <deal.II/base/multithread_info.h>
21 #include <deal.II/base/parallel.h>
22 #include <deal.II/base/utilities.h>
23 
24 #include <deal.II/matrix_free/task_info.h>
25 
26 
27 #ifdef DEAL_II_WITH_THREADS
28 # include <tbb/blocked_range.h>
29 # include <tbb/parallel_for.h>
30 # include <tbb/task.h>
31 # include <tbb/task_scheduler_init.h>
32 #endif
33 
34 #include <iostream>
35 #include <set>
36 
37 DEAL_II_NAMESPACE_OPEN
38 
39 
40 
41 /*-------------------- Implementation of the matrix-free loop --------------*/
42 namespace internal
43 {
44  namespace MatrixFreeFunctions
45  {
46 #ifdef DEAL_II_WITH_THREADS
47 
48  // This defines the TBB data structures that are needed to schedule the
49  // partition-partition variant
50 
51  namespace partition
52  {
53  class ActualCellWork
54  {
55  public:
56  ActualCellWork(MFWorkerInterface **worker_pointer,
57  const unsigned int partition,
58  const TaskInfo & task_info)
59  : worker(nullptr)
60  , worker_pointer(worker_pointer)
62  , task_info(task_info)
63  {}
64 
65  ActualCellWork(MFWorkerInterface &worker,
66  const unsigned int partition,
67  const TaskInfo & task_info)
68  : worker(&worker)
69  , worker_pointer(nullptr)
71  , task_info(task_info)
72  {}
73 
74  void
75  operator()() const
76  {
77  MFWorkerInterface *used_worker =
78  worker != nullptr ? worker : *worker_pointer;
79  Assert(used_worker != nullptr, ExcInternalError());
80  used_worker->cell(
81  std::make_pair(task_info.cell_partition_data[partition],
82  task_info.cell_partition_data[partition + 1]));
83 
84  if (task_info.face_partition_data.empty() == false)
85  {
86  used_worker->face(
87  std::make_pair(task_info.face_partition_data[partition],
88  task_info.face_partition_data[partition + 1]));
89 
90  used_worker->boundary(std::make_pair(
91  task_info.boundary_partition_data[partition],
92  task_info.boundary_partition_data[partition + 1]));
93  }
94  }
95 
96  private:
97  MFWorkerInterface * worker;
98  MFWorkerInterface **worker_pointer;
99  const unsigned int partition;
100  const TaskInfo & task_info;
101  };
102 
103  class CellWork : public tbb::task
104  {
105  public:
106  CellWork(MFWorkerInterface &worker,
107  const unsigned int partition,
108  const TaskInfo & task_info,
109  const bool is_blocked)
110  : dummy(nullptr)
111  , work(worker, partition, task_info)
112  , is_blocked(is_blocked)
113  {}
114 
115  tbb::task *
116  execute() override
117  {
118  work();
119 
120  if (is_blocked == true)
121  tbb::empty_task::spawn(*dummy);
122  return nullptr;
123  }
124 
125  tbb::empty_task *dummy;
126 
127  private:
128  ActualCellWork work;
129  const bool is_blocked;
130  };
131 
132 
133 
134  class PartitionWork : public tbb::task
135  {
136  public:
137  PartitionWork(MFWorkerInterface &function_in,
138  const unsigned int partition_in,
139  const TaskInfo & task_info_in,
140  const bool is_blocked_in = false)
141  : dummy(nullptr)
142  , function(function_in)
143  , partition(partition_in)
144  , task_info(task_info_in)
145  , is_blocked(is_blocked_in)
146  {}
147 
148  tbb::task *
149  execute() override
150  {
151  tbb::empty_task *root =
152  new (tbb::task::allocate_root()) tbb::empty_task;
153  const unsigned int evens = task_info.partition_evens[partition];
154  const unsigned int odds = task_info.partition_odds[partition];
155  const unsigned int n_blocked_workers =
156  task_info.partition_n_blocked_workers[partition];
157  const unsigned int n_workers =
158  task_info.partition_n_workers[partition];
159  std::vector<CellWork *> worker(n_workers);
160  std::vector<CellWork *> blocked_worker(n_blocked_workers);
161 
162  root->set_ref_count(evens + 1);
163  for (unsigned int j = 0; j < evens; j++)
164  {
165  worker[j] = new (root->allocate_child())
166  CellWork(function,
167  task_info.partition_row_index[partition] + 2 * j,
168  task_info,
169  false);
170  if (j > 0)
171  {
172  worker[j]->set_ref_count(2);
173  blocked_worker[j - 1]->dummy =
174  new (worker[j]->allocate_child()) tbb::empty_task;
175  tbb::task::spawn(*blocked_worker[j - 1]);
176  }
177  else
178  worker[j]->set_ref_count(1);
179  if (j < evens - 1)
180  {
181  blocked_worker[j] = new (worker[j]->allocate_child())
182  CellWork(function,
183  task_info.partition_row_index[partition] + 2 * j +
184  1,
185  task_info,
186  true);
187  }
188  else
189  {
190  if (odds == evens)
191  {
192  worker[evens] = new (worker[j]->allocate_child())
193  CellWork(function,
194  task_info.partition_row_index[partition] +
195  2 * j + 1,
196  task_info,
197  false);
198  tbb::task::spawn(*worker[evens]);
199  }
200  else
201  {
202  tbb::empty_task *child =
203  new (worker[j]->allocate_child()) tbb::empty_task();
204  tbb::task::spawn(*child);
205  }
206  }
207  }
208 
209  root->wait_for_all();
210  root->destroy(*root);
211  if (is_blocked == true)
212  tbb::empty_task::spawn(*dummy);
213  return nullptr;
214  }
215 
216  tbb::empty_task *dummy;
217 
218  private:
219  MFWorkerInterface &function;
220  const unsigned int partition;
221  const TaskInfo & task_info;
222  const bool is_blocked;
223  };
224 
225  } // end of namespace partition
226 
227 
228 
229  namespace color
230  {
231  class CellWork
232  {
233  public:
234  CellWork(MFWorkerInterface &worker_in,
235  const TaskInfo & task_info_in,
236  const unsigned int partition_in)
237  : worker(worker_in)
238  , task_info(task_info_in)
239  , partition(partition_in)
240  {}
241 
242  void
243  operator()(const tbb::blocked_range<unsigned int> &r) const
244  {
245  const unsigned int start_index =
246  task_info.cell_partition_data[partition] +
247  task_info.block_size * r.begin();
248  const unsigned int end_index =
249  std::min(start_index + task_info.block_size * (r.end() - r.begin()),
250  task_info.cell_partition_data[partition + 1]);
251  worker.cell(std::make_pair(start_index, end_index));
252 
253  if (task_info.face_partition_data.empty() == false)
254  {
255  AssertThrow(false, ExcNotImplemented());
256  }
257  }
258 
259  private:
260  MFWorkerInterface &worker;
261  const TaskInfo & task_info;
262  const unsigned int partition;
263  };
264 
265 
266 
267  class PartitionWork : public tbb::task
268  {
269  public:
270  PartitionWork(MFWorkerInterface &worker_in,
271  const unsigned int partition_in,
272  const TaskInfo & task_info_in,
273  const bool is_blocked_in)
274  : dummy(nullptr)
275  , worker(worker_in)
276  , partition(partition_in)
277  , task_info(task_info_in)
278  , is_blocked(is_blocked_in)
279  {}
280 
281  tbb::task *
282  execute() override
283  {
284  const unsigned int n_chunks =
285  (task_info.cell_partition_data[partition + 1] -
286  task_info.cell_partition_data[partition] + task_info.block_size -
287  1) /
288  task_info.block_size;
289  parallel_for(tbb::blocked_range<unsigned int>(0, n_chunks, 1),
290  CellWork(worker, task_info, partition));
291  if (is_blocked == true)
292  tbb::empty_task::spawn(*dummy);
293  return nullptr;
294  }
295 
296  tbb::empty_task *dummy;
297 
298  private:
299  MFWorkerInterface &worker;
300  const unsigned int partition;
301  const TaskInfo & task_info;
302  const bool is_blocked;
303  };
304 
305  } // end of namespace color
306 
307 
308 
309  class MPICommunication : public tbb::task
310  {
311  public:
312  MPICommunication(MFWorkerInterface &worker_in, const bool do_compress)
313  : worker(worker_in)
314  , do_compress(do_compress)
315  {}
316 
317  tbb::task *
318  execute() override
319  {
320  if (do_compress == false)
321  worker.vector_update_ghosts_finish();
322  else
323  worker.vector_compress_start();
324  return nullptr;
325  }
326 
327  private:
328  MFWorkerInterface &worker;
329  const bool do_compress;
330  };
331 
332 #endif // DEAL_II_WITH_THREADS
333 
334 
335 
336  void
338  {
340 
341 #ifdef DEAL_II_WITH_THREADS
342 
343  if (scheme != none)
344  {
346  if (scheme == partition_partition)
347  {
348  tbb::empty_task *root =
349  new (tbb::task::allocate_root()) tbb::empty_task;
350  root->set_ref_count(evens + 1);
351  std::vector<partition::PartitionWork *> worker(n_workers);
352  std::vector<partition::PartitionWork *> blocked_worker(
354  MPICommunication *worker_compr =
355  new (root->allocate_child()) MPICommunication(funct, true);
356  worker_compr->set_ref_count(1);
357  for (unsigned int j = 0; j < evens; j++)
358  {
359  if (j > 0)
360  {
361  worker[j] = new (root->allocate_child())
362  partition::PartitionWork(funct, 2 * j, *this, false);
363  worker[j]->set_ref_count(2);
364  blocked_worker[j - 1]->dummy =
365  new (worker[j]->allocate_child()) tbb::empty_task;
366  tbb::task::spawn(*blocked_worker[j - 1]);
367  }
368  else
369  {
370  worker[j] = new (worker_compr->allocate_child())
371  partition::PartitionWork(funct, 2 * j, *this, false);
372  worker[j]->set_ref_count(2);
373  MPICommunication *worker_dist =
374  new (worker[j]->allocate_child())
375  MPICommunication(funct, false);
376  tbb::task::spawn(*worker_dist);
377  }
378  if (j < evens - 1)
379  {
380  blocked_worker[j] = new (worker[j]->allocate_child())
381  partition::PartitionWork(funct, 2 * j + 1, *this, true);
382  }
383  else
384  {
385  if (odds == evens)
386  {
387  worker[evens] = new (worker[j]->allocate_child())
388  partition::PartitionWork(funct,
389  2 * j + 1,
390  *this,
391  false);
392  tbb::task::spawn(*worker[evens]);
393  }
394  else
395  {
396  tbb::empty_task *child =
397  new (worker[j]->allocate_child()) tbb::empty_task();
398  tbb::task::spawn(*child);
399  }
400  }
401  }
402 
403  root->wait_for_all();
404  root->destroy(*root);
405  }
406  else // end of partition-partition, start of partition-color
407  {
408  // check whether there is only one partition. if not, build up the
409  // tree of partitions
410  if (odds > 0)
411  {
412  tbb::empty_task *root =
413  new (tbb::task::allocate_root()) tbb::empty_task;
414  root->set_ref_count(evens + 1);
415  const unsigned int n_blocked_workers =
416  odds - (odds + evens + 1) % 2;
417  const unsigned int n_workers =
419  std::vector<color::PartitionWork *> worker(n_workers);
420  std::vector<color::PartitionWork *> blocked_worker(
422  unsigned int worker_index = 0, slice_index = 0;
423  int spawn_index_child = -2;
424  MPICommunication *worker_compr =
425  new (root->allocate_child()) MPICommunication(funct, true);
426  worker_compr->set_ref_count(1);
427  for (unsigned int part = 0;
428  part < partition_row_index.size() - 1;
429  part++)
430  {
431  if (part == 0)
432  worker[worker_index] =
433  new (worker_compr->allocate_child())
434  color::PartitionWork(funct,
435  slice_index,
436  *this,
437  false);
438  else
439  worker[worker_index] = new (root->allocate_child())
440  color::PartitionWork(funct,
441  slice_index,
442  *this,
443  false);
444  slice_index++;
445  for (; slice_index < partition_row_index[part + 1];
446  slice_index++)
447  {
448  worker[worker_index]->set_ref_count(1);
449  worker_index++;
450  worker[worker_index] =
451  new (worker[worker_index - 1]->allocate_child())
452  color::PartitionWork(funct,
453  slice_index,
454  *this,
455  false);
456  }
457  worker[worker_index]->set_ref_count(2);
458  if (part > 0)
459  {
460  blocked_worker[(part - 1) / 2]->dummy =
461  new (worker[worker_index]->allocate_child())
462  tbb::empty_task;
463  worker_index++;
464  if (spawn_index_child == -1)
465  tbb::task::spawn(*blocked_worker[(part - 1) / 2]);
466  else
467  {
468  Assert(spawn_index_child >= 0,
469  ExcInternalError());
470  tbb::task::spawn(*worker[spawn_index_child]);
471  }
472  spawn_index_child = -2;
473  }
474  else
475  {
476  MPICommunication *worker_dist =
477  new (worker[worker_index]->allocate_child())
478  MPICommunication(funct, false);
479  tbb::task::spawn(*worker_dist);
480  worker_index++;
481  }
482  part += 1;
483  if (part < partition_row_index.size() - 1)
484  {
485  if (part < partition_row_index.size() - 2)
486  {
487  blocked_worker[part / 2] =
488  new (worker[worker_index - 1]->allocate_child())
489  color::PartitionWork(funct,
490  slice_index,
491  *this,
492  true);
493  slice_index++;
494  if (slice_index < partition_row_index[part + 1])
495  {
496  blocked_worker[part / 2]->set_ref_count(1);
497  worker[worker_index] = new (
498  blocked_worker[part / 2]->allocate_child())
499  color::PartitionWork(funct,
500  slice_index,
501  *this,
502  false);
503  slice_index++;
504  }
505  else
506  {
507  spawn_index_child = -1;
508  continue;
509  }
510  }
511  for (; slice_index < partition_row_index[part + 1];
512  slice_index++)
513  {
514  if (slice_index > partition_row_index[part])
515  {
516  worker[worker_index]->set_ref_count(1);
517  worker_index++;
518  }
519  worker[worker_index] =
520  new (worker[worker_index - 1]->allocate_child())
521  color::PartitionWork(funct,
522  slice_index,
523  *this,
524  false);
525  }
526  spawn_index_child = worker_index;
527  worker_index++;
528  }
529  else
530  {
531  tbb::empty_task *final =
532  new (worker[worker_index - 1]->allocate_child())
533  tbb::empty_task;
534  tbb::task::spawn(*final);
535  spawn_index_child = worker_index - 1;
536  }
537  }
538  if (evens == odds)
539  {
540  Assert(spawn_index_child >= 0, ExcInternalError());
541  tbb::task::spawn(*worker[spawn_index_child]);
542  }
543  root->wait_for_all();
544  root->destroy(*root);
545  }
546  // case when we only have one partition: this is the usual
547  // coloring scheme, and we just schedule a parallel for loop for
548  // each color
549  else
550  {
551  Assert(evens <= 1, ExcInternalError());
553 
554  for (unsigned int color = 0; color < partition_row_index[1];
555  ++color)
556  {
557  tbb::empty_task *root =
558  new (tbb::task::allocate_root()) tbb::empty_task;
559  root->set_ref_count(2);
560  color::PartitionWork *worker =
561  new (root->allocate_child())
562  color::PartitionWork(funct, color, *this, false);
563  tbb::empty_task::spawn(*worker);
564  root->wait_for_all();
565  root->destroy(*root);
566  }
567 
568  funct.vector_compress_start();
569  }
570  }
571  }
572  else
573 #endif
574  // serial loop, go through up to three times and do the MPI transfer at
575  // the beginning/end of the second part
576  {
577  for (unsigned int part = 0; part < partition_row_index.size() - 2;
578  ++part)
579  {
580  if (part == 1)
582 
583  for (unsigned int i = partition_row_index[part];
584  i < partition_row_index[part + 1];
585  ++i)
586  {
587  AssertIndexRange(i + 1, cell_partition_data.size());
589  {
590  funct.zero_dst_vector_range(i);
591  funct.cell(std::make_pair(cell_partition_data[i],
592  cell_partition_data[i + 1]));
593  }
594 
595  if (face_partition_data.empty() == false)
596  {
598  funct.face(std::make_pair(face_partition_data[i],
599  face_partition_data[i + 1]));
600  if (boundary_partition_data[i + 1] >
602  funct.boundary(
603  std::make_pair(boundary_partition_data[i],
604  boundary_partition_data[i + 1]));
605  }
606  }
607 
608  if (part == 1)
609  funct.vector_compress_start();
610  }
611  }
612  funct.vector_compress_finish();
613  }
614 
615 
616 
618  {
619  clear();
620  }
621 
622 
623 
624  void
626  {
627  n_active_cells = 0;
628  n_ghost_cells = 0;
630  block_size = 0;
631  n_blocks = 0;
632  scheme = none;
633  partition_row_index.clear();
634  partition_row_index.resize(2);
635  cell_partition_data.clear();
636  face_partition_data.clear();
637  boundary_partition_data.clear();
638  evens = 0;
639  odds = 0;
640  n_blocked_workers = 0;
641  n_workers = 0;
642  partition_evens.clear();
643  partition_odds.clear();
645  partition_n_workers.clear();
646  communicator = MPI_COMM_SELF;
647  my_pid = 0;
648  n_procs = 1;
649  }
650 
651 
652 
653  template <typename StreamType>
654  void
656  const std::size_t data_length) const
657  {
658  Utilities::MPI::MinMaxAvg memory_c =
659  Utilities::MPI::min_max_avg(1e-6 * data_length, communicator);
660  if (n_procs < 2)
661  out << memory_c.min;
662  else
663  out << memory_c.min << "/" << memory_c.avg << "/" << memory_c.max;
664  out << " MB" << std::endl;
665  }
666 
667 
668 
669  std::size_t
671  {
672  return (
673  sizeof(*this) +
682  }
683 
684 
685 
686  void
688  const unsigned int n_active_cells_in,
689  const unsigned int n_active_and_ghost_cells,
690  const unsigned int vectorization_length_in,
691  std::vector<unsigned int> &boundary_cells)
692  {
693  vectorization_length = vectorization_length_in;
694  n_active_cells = n_active_cells_in;
695  n_ghost_cells = n_active_and_ghost_cells - n_active_cells;
696 
697  // try to make the number of boundary cells divisible by the number of
698  // vectors in vectorization
699  unsigned int fillup_needed =
700  (vectorization_length - boundary_cells.size() % vectorization_length) %
702  if (fillup_needed > 0 && boundary_cells.size() < n_active_cells)
703  {
704  // fill additional cells into the list of boundary cells to get a
705  // balanced number. Go through the indices successively until we
706  // found enough indices
707  std::vector<unsigned int> new_boundary_cells;
708  new_boundary_cells.reserve(boundary_cells.size());
709 
710  unsigned int next_free_slot = 0, bound_index = 0;
711  while (fillup_needed > 0 && bound_index < boundary_cells.size())
712  {
713  if (next_free_slot < boundary_cells[bound_index])
714  {
715  // check if there are enough cells to fill with in the
716  // current slot
717  if (next_free_slot + fillup_needed <=
718  boundary_cells[bound_index])
719  {
720  for (unsigned int j =
721  boundary_cells[bound_index] - fillup_needed;
722  j < boundary_cells[bound_index];
723  ++j)
724  new_boundary_cells.push_back(j);
725  fillup_needed = 0;
726  }
727  // ok, not enough indices, so just take them all up to the
728  // next boundary cell
729  else
730  {
731  for (unsigned int j = next_free_slot;
732  j < boundary_cells[bound_index];
733  ++j)
734  new_boundary_cells.push_back(j);
735  fillup_needed -=
736  boundary_cells[bound_index] - next_free_slot;
737  }
738  }
739  new_boundary_cells.push_back(boundary_cells[bound_index]);
740  next_free_slot = boundary_cells[bound_index] + 1;
741  ++bound_index;
742  }
743  while (fillup_needed > 0 &&
744  (new_boundary_cells.size() == 0 ||
745  new_boundary_cells.back() < n_active_cells - 1))
746  new_boundary_cells.push_back(new_boundary_cells.back() + 1);
747  while (bound_index < boundary_cells.size())
748  new_boundary_cells.push_back(boundary_cells[bound_index++]);
749 
750  boundary_cells.swap(new_boundary_cells);
751  }
752 
753  // set the number of cells
754  std::sort(boundary_cells.begin(), boundary_cells.end());
755 
756  // check that number of boundary cells is divisible by
757  // vectorization_length or that it contains all cells
758  Assert(boundary_cells.size() % vectorization_length == 0 ||
759  boundary_cells.size() == n_active_cells,
760  ExcInternalError());
761  }
762 
763 
764 
765  void
767  const std::vector<unsigned int> &boundary_cells,
768  const std::vector<unsigned int> &cells_close_to_boundary,
769  const unsigned int dofs_per_cell,
770  const std::vector<unsigned int> &cell_vectorization_categories,
771  const bool cell_vectorization_categories_strict,
772  std::vector<unsigned int> & renumbering,
773  std::vector<unsigned char> & incompletely_filled_vectorization)
774  {
775  const unsigned int n_macro_cells =
777  const unsigned int n_ghost_slots =
779  const unsigned int n_boundary_cells = boundary_cells.size();
780 
781  incompletely_filled_vectorization.resize(n_macro_cells + n_ghost_slots);
782  renumbering.resize(n_active_cells + n_ghost_cells,
784 
785  // Define the outer number of partitions. In the MPI case, we have three
786  // partitions (part before comm, part with comm, part after comm)
787  if (n_procs == 1)
788  partition_row_index.resize(3);
789  else
790  partition_row_index.resize(5);
791 
792  // Initially mark the cells according to the MPI ranking
793  std::vector<unsigned char> cell_marked(n_active_cells + n_ghost_cells, 0);
794  if (n_procs > 1)
795  {
796  for (unsigned int i = 0; i < n_boundary_cells; ++i)
797  cell_marked[boundary_cells[i]] = 2;
798 
799  Assert(boundary_cells.size() % vectorization_length == 0 ||
800  boundary_cells.size() == n_active_cells,
801  ExcInternalError());
802 
803  const unsigned int n_second_slot =
804  ((n_active_cells - n_boundary_cells) / 2 / vectorization_length) *
806  unsigned int count = 0;
807  for (const unsigned int cell : cells_close_to_boundary)
808  if (cell_marked[cell] == 0)
809  {
810  cell_marked[cell] = count < n_second_slot ? 1 : 3;
811  ++count;
812  }
813 
814  unsigned int c = 0;
815  for (; c < n_active_cells && count < n_second_slot; ++c)
816  if (cell_marked[c] == 0)
817  {
818  cell_marked[c] = 1;
819  ++count;
820  }
821  for (; c < n_active_cells; ++c)
822  if (cell_marked[c] == 0)
823  cell_marked[c] = 3;
824  for (; c < n_active_cells + n_ghost_cells; ++c)
825  if (cell_marked[c] == 0)
826  cell_marked[c] = 4;
827  }
828  else
829  std::fill(cell_marked.begin(), cell_marked.end(), 1);
830 
831  for (const unsigned char marker : cell_marked)
832  {
833  (void)marker;
834  Assert(marker != 0, ExcInternalError());
835  }
836 
837  unsigned int n_categories = 1;
838  std::vector<unsigned int> tight_category_map;
839  if (cell_vectorization_categories.empty() == false)
840  {
841  AssertDimension(cell_vectorization_categories.size(),
843 
844  // create a tight map of categories for not taking exceeding amounts
845  // of memory below. Sort the new categories by the numbers in the
846  // old one.
847  tight_category_map.resize(n_active_cells + n_ghost_cells);
848  std::set<unsigned int> used_categories;
849  for (unsigned int i = 0; i < n_active_cells + n_ghost_cells; ++i)
850  used_categories.insert(cell_vectorization_categories[i]);
851  std::vector<unsigned int> used_categories_vector(
852  used_categories.size());
853  n_categories = 0;
854  for (const auto &it : used_categories)
855  used_categories_vector[n_categories++] = it;
856  for (unsigned int i = 0; i < n_active_cells + n_ghost_cells; ++i)
857  {
858  const unsigned int index =
859  std::lower_bound(used_categories_vector.begin(),
860  used_categories_vector.end(),
861  cell_vectorization_categories[i]) -
862  used_categories_vector.begin();
863  AssertIndexRange(index, used_categories_vector.size());
864  tight_category_map[i] = index;
865  }
866 
867  // leave some more space for empty lanes
868  incompletely_filled_vectorization.resize(
869  incompletely_filled_vectorization.size() + 4 * n_categories);
870  }
871  else if (cells_close_to_boundary.empty())
872  tight_category_map.resize(n_active_cells + n_ghost_cells, 0);
873  else
874  {
875  n_categories = 2;
876  tight_category_map.resize(n_active_cells + n_ghost_cells, 1);
877  for (const unsigned int cell : cells_close_to_boundary)
878  tight_category_map[cell] = 0;
879  }
880 
881  cell_partition_data.clear();
882  cell_partition_data.resize(1, 0);
883  unsigned int counter = 0;
884  unsigned int n_cells = 0;
885  std::vector<std::vector<unsigned int>> renumbering_category(n_categories);
886  for (unsigned int block = 1; block < (n_procs > 1u ? 5u : 3u); ++block)
887  {
888  // step 1: sort by category
889  for (unsigned int i = 0; i < n_active_cells + n_ghost_cells; ++i)
890  if (cell_marked[i] == block)
891  renumbering_category[tight_category_map[i]].push_back(i);
892 
893  // step 2: if we want to fill up the ranges in vectorization, promote
894  // some of the cells to a higher category
895  if (cell_vectorization_categories_strict == false && n_categories > 1)
896  for (unsigned int j = n_categories - 1; j > 0; --j)
897  {
898  unsigned int lower_index = j - 1;
899  while (renumbering_category[j].size() % vectorization_length)
900  {
901  while (renumbering_category[j].size() %
903  !renumbering_category[lower_index].empty())
904  {
905  renumbering_category[j].push_back(
906  renumbering_category[lower_index].back());
907  renumbering_category[lower_index].pop_back();
908  }
909  if (lower_index == 0)
910  break;
911  else
912  --lower_index;
913  }
914  }
915 
916  // step 3: append cells according to categories
917  for (unsigned int j = 0; j < n_categories; ++j)
918  {
919  for (const unsigned int cell : renumbering_category[j])
920  renumbering[counter++] = cell;
921  unsigned int remainder =
922  renumbering_category[j].size() % vectorization_length;
923  if (remainder)
924  incompletely_filled_vectorization
925  [renumbering_category[j].size() / vectorization_length +
926  n_cells] = remainder;
927  const unsigned int n_my_macro_cells =
928  (renumbering_category[j].size() + vectorization_length - 1) /
930  renumbering_category[j].clear();
931 
932  // step 4: create blocks for face integrals, make the number of
933  // cells divisible by 4 if possible
934  const unsigned int block_size =
935  std::max((2048U / dofs_per_cell) / 8 * 4, 2U);
936  if (block < 4)
937  for (unsigned int k = 0; k < n_my_macro_cells; k += block_size)
938  cell_partition_data.push_back(
939  n_cells + std::min(k + block_size, n_my_macro_cells));
940  else
941  cell_partition_data.back() += n_my_macro_cells;
942  n_cells += n_my_macro_cells;
943  }
944  partition_row_index[block] = cell_partition_data.size() - 1;
945  if (block == 3 || (block == 1 && n_procs == 1))
946  cell_partition_data.push_back(n_cells);
947  }
948  if (cell_vectorization_categories_strict == true)
949  {
950  Assert(n_cells >= n_macro_cells + n_ghost_slots, ExcInternalError());
951  }
952  else
953  {
954  AssertDimension(n_cells, n_macro_cells + n_ghost_slots);
955  }
956  AssertDimension(cell_partition_data.back(), n_cells);
958 
959  incompletely_filled_vectorization.resize(cell_partition_data.back());
960  }
961 
962 
963 
964  void
966  const std::vector<unsigned int> &boundary_cells,
967  std::vector<unsigned int> & renumbering,
968  std::vector<unsigned char> & incompletely_filled_vectorization)
969  {
970  const unsigned int n_macro_cells =
972  const unsigned int n_ghost_slots =
974  incompletely_filled_vectorization.resize(n_macro_cells + n_ghost_slots);
975  if (n_macro_cells * vectorization_length > n_active_cells)
976  incompletely_filled_vectorization[n_macro_cells - 1] =
978  (n_macro_cells * vectorization_length - n_active_cells);
979  if (n_ghost_slots * vectorization_length > n_ghost_cells)
980  incompletely_filled_vectorization[n_macro_cells + n_ghost_slots - 1] =
982  (n_ghost_slots * vectorization_length - n_ghost_cells);
983 
984  std::vector<unsigned int> reverse_numbering(
986  for (unsigned int j = 0; j < boundary_cells.size(); ++j)
987  reverse_numbering[boundary_cells[j]] = j;
988  unsigned int counter = boundary_cells.size();
989  for (unsigned int j = 0; j < n_active_cells; ++j)
990  if (reverse_numbering[j] == numbers::invalid_unsigned_int)
991  reverse_numbering[j] = counter++;
992 
994  renumbering = Utilities::invert_permutation(reverse_numbering);
995 
996  for (unsigned int j = n_active_cells; j < n_active_cells + n_ghost_cells;
997  ++j)
998  renumbering.push_back(j);
999 
1000  // TODO: might be able to simplify this code by not relying on the cell
1001  // partition data while computing the thread graph
1002  cell_partition_data.clear();
1003  cell_partition_data.push_back(0);
1004  if (n_procs > 1)
1005  {
1006  const unsigned int n_macro_boundary_cells =
1007  (boundary_cells.size() + vectorization_length - 1) /
1009  cell_partition_data.push_back(
1010  (n_macro_cells - n_macro_boundary_cells) / 2);
1012  n_macro_boundary_cells);
1013  }
1014  else
1015  AssertDimension(boundary_cells.size(), 0);
1016  cell_partition_data.push_back(n_macro_cells);
1017  cell_partition_data.push_back(cell_partition_data.back() + n_ghost_slots);
1018  partition_row_index.resize(n_procs > 1 ? 4 : 2);
1019  partition_row_index[0] = 0;
1020  partition_row_index[1] = 1;
1021  if (n_procs > 1)
1022  {
1023  partition_row_index[2] = 2;
1024  partition_row_index[3] = 3;
1025  }
1026  }
1027 
1028 
1029 
1030  void
1031  TaskInfo::guess_block_size(const unsigned int dofs_per_cell)
1032  {
1033  // user did not say a positive number, so we have to guess
1034  if (block_size == 0)
1035  {
1036  // we would like to have enough work to do, so as first guess, try
1037  // to get 16 times as many chunks as we have threads on the system.
1040 
1041  // if there are too few degrees of freedom per cell, need to
1042  // increase the block size
1043  const unsigned int minimum_parallel_grain_size = 200;
1044  if (dofs_per_cell * block_size < minimum_parallel_grain_size)
1045  block_size = (minimum_parallel_grain_size / dofs_per_cell + 1);
1046  if (dofs_per_cell * block_size > 10000)
1047  block_size /= 4;
1048 
1049  block_size =
1050  1 << static_cast<unsigned int>(std::log2(block_size + 1));
1051  }
1052  if (block_size > n_active_cells)
1053  block_size = std::max(1U, n_active_cells);
1054  }
1055 
1056 
1057 
1058  void
1060  DynamicSparsityPattern & connectivity_large,
1061  std::vector<unsigned int> & renumbering,
1062  std::vector<unsigned char> &irregular_cells,
1063  const bool)
1064  {
1065  const unsigned int n_macro_cells = *(cell_partition_data.end() - 2);
1066  if (n_macro_cells == 0)
1067  return;
1068 
1070 
1071  unsigned int partition = 0, counter = 0;
1072 
1073  // Create connectivity graph for blocks based on connectivity graph for
1074  // cells.
1075  DynamicSparsityPattern connectivity(n_blocks, n_blocks);
1076  make_connectivity_cells_to_blocks(irregular_cells,
1077  connectivity_large,
1078  connectivity);
1079 
1080  // Create cell-block partitioning.
1081 
1082  // For each block of cells, this variable saves to which partitions the
1083  // block belongs. Initialize all to -1 to mark them as not yet assigned
1084  // a partition.
1085  std::vector<unsigned int> cell_partition(n_blocks,
1087 
1088  // In element j of this variable, one puts the old number of the block
1089  // that should be the jth block in the new numeration.
1090  std::vector<unsigned int> partition_list(n_blocks, 0);
1091  std::vector<unsigned int> partition_color_list(n_blocks, 0);
1092 
1093  // This vector points to the start of each partition.
1094  std::vector<unsigned int> partition_size(2, 0);
1095 
1096  // blocking_connectivity = true;
1097 
1098  // The cluster_size in make_partitioning defines that the no. of cells
1099  // in each partition should be a multiple of cluster_size.
1100  unsigned int cluster_size = 1;
1101 
1102  // Make the partitioning of the first layer of the blocks of cells.
1103  make_partitioning(connectivity,
1104  cluster_size,
1105  cell_partition,
1106  partition_list,
1107  partition_size,
1108  partition);
1109 
1110  // Color the cells within each partition
1112  partition,
1113  cell_partition,
1114  partition_list,
1115  partition_size,
1116  partition_color_list);
1117 
1118  partition_list = renumbering;
1119 
1120 #ifdef DEBUG
1121  // in debug mode, check that the partition color list is one-to-one
1122  {
1123  std::vector<unsigned int> sorted_pc_list(partition_color_list);
1124  std::sort(sorted_pc_list.begin(), sorted_pc_list.end());
1125  for (unsigned int i = 0; i < sorted_pc_list.size(); ++i)
1126  Assert(sorted_pc_list[i] == i, ExcInternalError());
1127  }
1128 #endif
1129 
1130  // set the start list for each block and compute the renumbering of
1131  // cells
1132  std::vector<unsigned int> block_start(n_macro_cells + 1);
1133  std::vector<unsigned char> irregular(n_macro_cells);
1134 
1135  unsigned int mcell_start = 0;
1136  block_start[0] = 0;
1137  for (unsigned int block = 0; block < n_blocks; block++)
1138  {
1139  block_start[block + 1] = block_start[block];
1140  for (unsigned int mcell = mcell_start;
1141  mcell < std::min(mcell_start + block_size, n_macro_cells);
1142  ++mcell)
1143  {
1144  unsigned int n_comp = (irregular_cells[mcell] > 0) ?
1145  irregular_cells[mcell] :
1147  block_start[block + 1] += n_comp;
1148  ++counter;
1149  }
1150  mcell_start += block_size;
1151  }
1152  counter = 0;
1153  unsigned int counter_macro = 0;
1154  unsigned int block_size_last =
1155  n_macro_cells - block_size * (n_blocks - 1);
1156  if (block_size_last == 0)
1157  block_size_last = block_size;
1158 
1159  unsigned int tick = 0;
1160  for (unsigned int block = 0; block < n_blocks; block++)
1161  {
1162  unsigned int present_block = partition_color_list[block];
1163  for (unsigned int cell = block_start[present_block];
1164  cell < block_start[present_block + 1];
1165  ++cell)
1166  renumbering[counter++] = partition_list[cell];
1167  unsigned int this_block_size =
1168  (present_block == n_blocks - 1) ? block_size_last : block_size;
1169 
1170  // Also re-compute the content of cell_partition_data to
1171  // contain the numbers of cells, not blocks
1172  if (cell_partition_data[tick] == block)
1173  cell_partition_data[tick++] = counter_macro;
1174 
1175  for (unsigned int j = 0; j < this_block_size; j++)
1176  irregular[counter_macro++] =
1177  irregular_cells[present_block * block_size + j];
1178  }
1179  AssertDimension(tick + 1, cell_partition_data.size());
1180  cell_partition_data.back() = counter_macro;
1181 
1182  irregular_cells.swap(irregular);
1183  AssertDimension(counter, n_active_cells);
1184  AssertDimension(counter_macro, n_macro_cells);
1185 
1186  // check that the renumbering is one-to-one
1187 #ifdef DEBUG
1188  {
1189  std::vector<unsigned int> sorted_renumbering(renumbering);
1190  std::sort(sorted_renumbering.begin(), sorted_renumbering.end());
1191  for (unsigned int i = 0; i < sorted_renumbering.size(); ++i)
1192  Assert(sorted_renumbering[i] == i, ExcInternalError());
1193  }
1194 #endif
1195 
1196 
1198  partition); // Actually sets too much for partition color case
1199 
1200  AssertDimension(cell_partition_data.back(), n_macro_cells);
1201  }
1202 
1203 
1204 
1205  void
1207  const std::vector<unsigned int> &cell_active_fe_index,
1208  DynamicSparsityPattern & connectivity,
1209  std::vector<unsigned int> & renumbering,
1210  std::vector<unsigned char> & irregular_cells,
1211  const bool hp_bool)
1212  {
1213  const unsigned int n_macro_cells = *(cell_partition_data.end() - 2);
1214  if (n_macro_cells == 0)
1215  return;
1216 
1218 
1219  // if we want to block before partitioning, create connectivity graph
1220  // for blocks based on connectivity graph for cells.
1221  DynamicSparsityPattern connectivity_blocks(n_blocks, n_blocks);
1222  make_connectivity_cells_to_blocks(irregular_cells,
1223  connectivity,
1224  connectivity_blocks);
1225 
1226  unsigned int n_blocks = 0;
1227  if (scheme == partition_color ||
1228  scheme == color) // blocking_connectivity == true
1229  n_blocks = this->n_blocks;
1230  else
1231  n_blocks = n_active_cells;
1232 
1233  // For each block of cells, this variable saves to which partitions the
1234  // block belongs. Initialize all to -1 to mark them as not yet assigned
1235  // a partition.
1236  std::vector<unsigned int> cell_partition(n_blocks,
1238 
1239  // In element j of this variable, one puts the old number (but after
1240  // renumbering according to the input renumbering) of the block that
1241  // should be the jth block in the new numeration.
1242  std::vector<unsigned int> partition_list(n_blocks, 0);
1243  std::vector<unsigned int> partition_2layers_list(n_blocks, 0);
1244 
1245  // This vector points to the start of each partition.
1246  std::vector<unsigned int> partition_size(2, 0);
1247 
1248  unsigned int partition = 0;
1249 
1250  // Within the partitions we want to be able to block for the case that
1251  // we do not block already in the connectivity. The cluster_size in
1252  // make_partitioning defines that the no. of cells in each partition
1253  // should be a multiple of cluster_size.
1254  unsigned int cluster_size = 1;
1255  if (scheme == partition_partition)
1256  cluster_size = block_size * vectorization_length;
1257 
1258  // Make the partitioning of the first layer of the blocks of cells.
1259  if (scheme == partition_color || scheme == color)
1260  make_partitioning(connectivity_blocks,
1261  cluster_size,
1262  cell_partition,
1263  partition_list,
1264  partition_size,
1265  partition);
1266  else
1267  make_partitioning(connectivity,
1268  cluster_size,
1269  cell_partition,
1270  partition_list,
1271  partition_size,
1272  partition);
1273 
1274  // Partition or color second layer
1275  if (scheme == partition_partition)
1276 
1277  {
1278  // Partition within partitions.
1280  connectivity,
1281  cell_active_fe_index,
1282  partition,
1283  cluster_size,
1284  hp_bool,
1285  cell_partition,
1286  partition_list,
1287  partition_size,
1288  partition_2layers_list,
1289  irregular_cells);
1290  }
1291  else if (scheme == partition_color || scheme == color)
1292  {
1293  make_coloring_within_partitions_pre_blocked(connectivity_blocks,
1294  partition,
1295  cell_partition,
1296  partition_list,
1297  partition_size,
1298  partition_2layers_list);
1299  }
1300 
1301  // in debug mode, check that the partition_2layers_list is one-to-one
1302 #ifdef DEBUG
1303  {
1304  std::vector<unsigned int> sorted_pc_list(partition_2layers_list);
1305  std::sort(sorted_pc_list.begin(), sorted_pc_list.end());
1306  for (unsigned int i = 0; i < sorted_pc_list.size(); ++i)
1307  Assert(sorted_pc_list[i] == i, ExcInternalError());
1308  }
1309 #endif
1310 
1311  // Set the new renumbering
1312  std::vector<unsigned int> renumbering_in(n_active_cells, 0);
1313  renumbering_in.swap(renumbering);
1314  if (scheme == partition_partition) // blocking_connectivity == false
1315  {
1316  // This is the simple case. The renumbering is just a combination of
1317  // the renumbering that we were given as an input and the
1318  // renumbering of partition/coloring given in partition_2layers_list
1319  for (unsigned int j = 0; j < renumbering.size(); j++)
1320  renumbering[j] = renumbering_in[partition_2layers_list[j]];
1321  // Account for the ghost cells, finally.
1322  for (unsigned int i = 0; i < n_ghost_cells; ++i)
1323  renumbering.push_back(i + n_active_cells);
1324  }
1325  else
1326  {
1327  // set the start list for each block and compute the renumbering of
1328  // cells
1329  std::vector<unsigned int> block_start(n_macro_cells + 1);
1330  std::vector<unsigned char> irregular(n_macro_cells);
1331 
1332  unsigned int counter = 0;
1333  unsigned int mcell_start = 0;
1334  block_start[0] = 0;
1335  for (unsigned int block = 0; block < n_blocks; block++)
1336  {
1337  block_start[block + 1] = block_start[block];
1338  for (unsigned int mcell = mcell_start;
1339  mcell < std::min(mcell_start + block_size, n_macro_cells);
1340  ++mcell)
1341  {
1342  unsigned int n_comp = (irregular_cells[mcell] > 0) ?
1343  irregular_cells[mcell] :
1345  block_start[block + 1] += n_comp;
1346  ++counter;
1347  }
1348  mcell_start += block_size;
1349  }
1350  counter = 0;
1351  unsigned int counter_macro = 0;
1352  unsigned int block_size_last =
1353  n_macro_cells - block_size * (n_blocks - 1);
1354  if (block_size_last == 0)
1355  block_size_last = block_size;
1356 
1357  unsigned int tick = 0;
1358  for (unsigned int block = 0; block < n_blocks; block++)
1359  {
1360  unsigned int present_block = partition_2layers_list[block];
1361  for (unsigned int cell = block_start[present_block];
1362  cell < block_start[present_block + 1];
1363  ++cell)
1364  renumbering[counter++] = renumbering_in[cell];
1365  unsigned int this_block_size =
1366  (present_block == n_blocks - 1) ? block_size_last : block_size;
1367 
1368  // Also re-compute the content of cell_partition_data to
1369  // contain the numbers of cells, not blocks
1370  if (cell_partition_data[tick] == block)
1371  cell_partition_data[tick++] = counter_macro;
1372 
1373  for (unsigned int j = 0; j < this_block_size; j++)
1374  irregular[counter_macro++] =
1375  irregular_cells[present_block * block_size + j];
1376  }
1377  AssertDimension(tick + 1, cell_partition_data.size());
1378  cell_partition_data.back() = counter_macro;
1379 
1380  irregular_cells.swap(irregular);
1381  AssertDimension(counter, n_active_cells);
1382  AssertDimension(counter_macro, n_macro_cells);
1383  // check that the renumbering is one-to-one
1384 #ifdef DEBUG
1385  {
1386  std::vector<unsigned int> sorted_renumbering(renumbering);
1387  std::sort(sorted_renumbering.begin(), sorted_renumbering.end());
1388  for (unsigned int i = 0; i < sorted_renumbering.size(); ++i)
1389  Assert(sorted_renumbering[i] == i, ExcInternalError());
1390  }
1391 #endif
1392  }
1393 
1394  // Update the task_info with the more information for the thread graph.
1395  update_task_info(partition);
1396  }
1397 
1398 
1399 
1400  void
1402  const std::vector<unsigned int> &cell_active_fe_index,
1403  DynamicSparsityPattern & connectivity,
1404  std::vector<unsigned int> & renumbering,
1405  std::vector<unsigned char> & irregular_cells,
1406  const bool hp_bool)
1407  {
1408  const unsigned int n_macro_cells = *(cell_partition_data.end() - 2);
1409  if (n_macro_cells == 0)
1410  return;
1411 
1412  const unsigned int cluster_size = block_size * vectorization_length;
1413 
1414  // Create cell-block partitioning.
1415 
1416  // For each block of cells, this variable saves to which partitions the
1417  // block belongs. Initialize all to n_macro_cells to mark them as not
1418  // yet assigned a partition.
1419  std::vector<unsigned int> cell_partition(n_active_cells,
1421 
1422 
1423  // In element j of this variable, one puts the old number of the block
1424  // that should be the jth block in the new numeration.
1425  std::vector<unsigned int> partition_list(n_active_cells, 0);
1426  std::vector<unsigned int> partition_partition_list(n_active_cells, 0);
1427 
1428  // This vector points to the start of each partition.
1429  std::vector<unsigned int> partition_size(2, 0);
1430 
1431  unsigned int partition = 0;
1432  // Here, we do not block inside the connectivity graph
1433  // blocking_connectivity = false;
1434 
1435  // Make the partitioning of the first layer of the blocks of cells.
1436  make_partitioning(connectivity,
1437  cluster_size,
1438  cell_partition,
1439  partition_list,
1440  partition_size,
1441  partition);
1442 
1443  // Partition within partitions.
1445  cell_active_fe_index,
1446  partition,
1447  cluster_size,
1448  hp_bool,
1449  cell_partition,
1450  partition_list,
1451  partition_size,
1452  partition_partition_list,
1453  irregular_cells);
1454 
1455  partition_list.swap(renumbering);
1456 
1457  for (unsigned int j = 0; j < renumbering.size(); j++)
1458  renumbering[j] = partition_list[partition_partition_list[j]];
1459 
1460  for (unsigned int i = 0; i < n_ghost_cells; ++i)
1461  renumbering.push_back(i + n_active_cells);
1462 
1463  update_task_info(partition);
1464  }
1465 
1466 
1467 
1468  void
1470  const std::vector<unsigned char> &irregular_cells,
1471  const DynamicSparsityPattern & connectivity_cells,
1472  DynamicSparsityPattern & connectivity_blocks) const
1473  {
1474  std::vector<std::vector<unsigned int>> cell_blocks(n_blocks);
1475  std::vector<unsigned int> touched_cells(n_active_cells);
1476  unsigned int cell = 0;
1477  for (unsigned int i = 0, mcell = 0; i < n_blocks; ++i)
1478  {
1479  for (unsigned int c = 0;
1480  c < block_size && mcell < *(cell_partition_data.end() - 2);
1481  ++c, ++mcell)
1482  {
1483  unsigned int ncomp = (irregular_cells[mcell] > 0) ?
1484  irregular_cells[mcell] :
1486  for (unsigned int c = 0; c < ncomp; ++c, ++cell)
1487  {
1488  cell_blocks[i].push_back(cell);
1489  touched_cells[cell] = i;
1490  }
1491  }
1492  }
1494  for (unsigned int i = 0; i < cell_blocks.size(); ++i)
1495  for (unsigned int col = 0; col < cell_blocks[i].size(); ++col)
1496  {
1498  connectivity_cells.begin(cell_blocks[i][col]);
1499  it != connectivity_cells.end(cell_blocks[i][col]);
1500  ++it)
1501  {
1502  if (touched_cells[it->column()] != i)
1503  connectivity_blocks.add(i, touched_cells[it->column()]);
1504  }
1505  }
1506  }
1507 
1508 
1509 
1510  // Function to create partitioning on the second layer within each
1511  // partition. Version without preblocking.
1512  void
1514  const DynamicSparsityPattern & connectivity,
1515  const std::vector<unsigned int> &cell_active_fe_index,
1516  const unsigned int partition,
1517  const unsigned int cluster_size,
1518  const bool hp_bool,
1519  const std::vector<unsigned int> &cell_partition,
1520  const std::vector<unsigned int> &partition_list,
1521  const std::vector<unsigned int> &partition_size,
1522  std::vector<unsigned int> & partition_partition_list,
1523  std::vector<unsigned char> & irregular_cells)
1524  {
1525  const unsigned int n_macro_cells = *(cell_partition_data.end() - 2);
1526  const unsigned int n_ghost_slots =
1527  *(cell_partition_data.end() - 1) - n_macro_cells;
1528 
1529  // List of cells in previous partition
1530  std::vector<unsigned int> neighbor_list;
1531  // List of cells in current partition for use as neighbors in next
1532  // partition
1533  std::vector<unsigned int> neighbor_neighbor_list;
1534 
1535  std::vector<unsigned int> renumbering(n_active_cells);
1536 
1537  irregular_cells.back() = 0;
1538  irregular_cells.resize(n_active_cells + n_ghost_slots);
1539 
1540  unsigned int max_fe_index = 0;
1541  for (const unsigned int fe_index : cell_active_fe_index)
1542  max_fe_index = std::max(fe_index, max_fe_index);
1543 
1544  Assert(!hp_bool || cell_active_fe_index.size() == n_active_cells,
1545  ExcInternalError());
1546 
1547  {
1548  unsigned int n_macro_cells_before = 0;
1549  // Create partitioning within partitions.
1550 
1551  // For each block of cells, this variable saves to which partitions
1552  // the block belongs. Initialize all to n_macro_cells to mark them as
1553  // not yet assigned a partition.
1554  std::vector<unsigned int> cell_partition_l2(
1556  partition_row_index.clear();
1557  partition_row_index.resize(partition + 1, 0);
1558  cell_partition_data.resize(1, 0);
1559 
1560  unsigned int counter = 0;
1561  unsigned int missing_macros;
1562  for (unsigned int part = 0; part < partition; ++part)
1563  {
1564  neighbor_neighbor_list.resize(0);
1565  neighbor_list.resize(0);
1566  bool work = true;
1567  unsigned int partition_l2 = 0;
1568  unsigned int start_up = partition_size[part];
1569  unsigned int partition_counter = 0;
1570  while (work)
1571  {
1572  if (neighbor_list.size() == 0)
1573  {
1574  work = false;
1575  partition_counter = 0;
1576  for (unsigned int j = start_up;
1577  j < partition_size[part + 1];
1578  ++j)
1579  if (cell_partition[partition_list[j]] == part &&
1580  cell_partition_l2[partition_list[j]] ==
1582  {
1583  start_up = j;
1584  work = true;
1585  partition_counter = 1;
1586  // To start up, set the start_up cell to partition
1587  // and list all its neighbors.
1588  AssertIndexRange(start_up, partition_size[part + 1]);
1589  cell_partition_l2[partition_list[start_up]] =
1590  partition_l2;
1591  neighbor_neighbor_list.push_back(
1592  partition_list[start_up]);
1593  partition_partition_list[counter++] =
1594  partition_list[start_up];
1595  start_up++;
1596  break;
1597  }
1598  }
1599  else
1600  {
1601  partition_counter = 0;
1602  for (const unsigned int neighbor : neighbor_list)
1603  {
1604  Assert(cell_partition[neighbor] == part,
1605  ExcInternalError());
1606  Assert(cell_partition_l2[neighbor] == partition_l2 - 1,
1607  ExcInternalError());
1608  auto neighbor_it = connectivity.begin(neighbor);
1609  const auto end_it = connectivity.end(neighbor);
1610  for (; neighbor_it != end_it; ++neighbor_it)
1611  {
1612  if (cell_partition[neighbor_it->column()] == part &&
1613  cell_partition_l2[neighbor_it->column()] ==
1615  {
1616  cell_partition_l2[neighbor_it->column()] =
1617  partition_l2;
1618  neighbor_neighbor_list.push_back(
1619  neighbor_it->column());
1620  partition_partition_list[counter++] =
1621  neighbor_it->column();
1622  partition_counter++;
1623  }
1624  }
1625  }
1626  }
1627  if (partition_counter > 0)
1628  {
1629  int index_before = neighbor_neighbor_list.size(),
1630  index = index_before;
1631  {
1632  // put the cells into separate lists for each FE index
1633  // within one partition-partition
1634  missing_macros = 0;
1635  std::vector<unsigned int> remaining_per_macro_cell(
1636  max_fe_index + 1);
1637  std::vector<std::vector<unsigned int>>
1638  renumbering_fe_index;
1639  unsigned int cell;
1640  bool filled = true;
1641  if (hp_bool == true)
1642  {
1643  renumbering_fe_index.resize(max_fe_index + 1);
1644  for (cell = counter - partition_counter;
1645  cell < counter;
1646  ++cell)
1647  {
1648  renumbering_fe_index
1649  [cell_active_fe_index.empty() ?
1650  0 :
1651  cell_active_fe_index
1652  [partition_partition_list[cell]]]
1653  .push_back(partition_partition_list[cell]);
1654  }
1655  // check how many more cells are needed in the lists
1656  for (unsigned int j = 0; j < max_fe_index + 1; j++)
1657  {
1658  remaining_per_macro_cell[j] =
1659  renumbering_fe_index[j].size() %
1661  if (remaining_per_macro_cell[j] != 0)
1662  filled = false;
1663  missing_macros +=
1664  ((renumbering_fe_index[j].size() +
1665  vectorization_length - 1) /
1667  }
1668  }
1669  else
1670  {
1671  remaining_per_macro_cell.resize(1);
1672  remaining_per_macro_cell[0] =
1673  partition_counter % vectorization_length;
1674  missing_macros =
1675  partition_counter / vectorization_length;
1676  if (remaining_per_macro_cell[0] != 0)
1677  {
1678  filled = false;
1679  missing_macros++;
1680  }
1681  }
1682  missing_macros =
1683  cluster_size - (missing_macros % cluster_size);
1684 
1685  // now we realized that there are some cells missing.
1686  while (missing_macros > 0 || filled == false)
1687  {
1688  if (index == 0)
1689  {
1690  index = neighbor_neighbor_list.size();
1691  if (index == index_before)
1692  {
1693  if (missing_macros != 0)
1694  {
1695  neighbor_neighbor_list.resize(0);
1696  }
1697  start_up--;
1698  break; // not connected - start again
1699  }
1700  index_before = index;
1701  }
1702  index--;
1703  unsigned int additional =
1704  neighbor_neighbor_list[index];
1705 
1706  // go through the neighbors of the last cell in the
1707  // current partition and check if we find some to
1708  // fill up with.
1710  connectivity.begin(
1711  additional),
1712  end =
1713  connectivity.end(
1714  additional);
1715  for (; neighbor != end; ++neighbor)
1716  {
1717  if (cell_partition[neighbor->column()] == part &&
1718  cell_partition_l2[neighbor->column()] ==
1720  {
1721  unsigned int this_index = 0;
1722  if (hp_bool == true)
1723  this_index =
1724  cell_active_fe_index.empty() ?
1725  0 :
1726  cell_active_fe_index[neighbor
1727  ->column()];
1728 
1729  // Only add this cell if we need more macro
1730  // cells in the current block or if there is
1731  // a macro cell with the FE index that is
1732  // not yet fully populated
1733  if (missing_macros > 0 ||
1734  remaining_per_macro_cell[this_index] > 0)
1735  {
1736  cell_partition_l2[neighbor->column()] =
1737  partition_l2;
1738  neighbor_neighbor_list.push_back(
1739  neighbor->column());
1740  if (hp_bool == true)
1741  renumbering_fe_index[this_index]
1742  .push_back(neighbor->column());
1743  partition_partition_list[counter] =
1744  neighbor->column();
1745  counter++;
1746  partition_counter++;
1747  if (remaining_per_macro_cell
1748  [this_index] == 0 &&
1749  missing_macros > 0)
1750  missing_macros--;
1751  remaining_per_macro_cell[this_index]++;
1752  if (remaining_per_macro_cell
1753  [this_index] ==
1755  {
1756  remaining_per_macro_cell[this_index] =
1757  0;
1758  }
1759  if (missing_macros == 0)
1760  {
1761  filled = true;
1762  for (unsigned int fe_ind = 0;
1763  fe_ind < max_fe_index + 1;
1764  ++fe_ind)
1765  if (remaining_per_macro_cell
1766  [fe_ind] != 0)
1767  filled = false;
1768  }
1769  if (filled == true)
1770  break;
1771  }
1772  }
1773  }
1774  }
1775  if (hp_bool == true)
1776  {
1777  // set the renumbering according to their active FE
1778  // index within one partition-partition which was
1779  // implicitly assumed above
1780  cell = counter - partition_counter;
1781  for (unsigned int j = 0; j < max_fe_index + 1; j++)
1782  {
1783  for (const unsigned int jj :
1784  renumbering_fe_index[j])
1785  renumbering[cell++] = jj;
1786  if (renumbering_fe_index[j].size() %
1788  0)
1789  irregular_cells[renumbering_fe_index[j].size() /
1791  n_macro_cells_before] =
1792  renumbering_fe_index[j].size() %
1794  n_macro_cells_before +=
1795  (renumbering_fe_index[j].size() +
1796  vectorization_length - 1) /
1798  renumbering_fe_index[j].resize(0);
1799  }
1800  }
1801  else
1802  {
1803  n_macro_cells_before +=
1804  partition_counter / vectorization_length;
1805  if (partition_counter % vectorization_length != 0)
1806  {
1807  irregular_cells[n_macro_cells_before] =
1808  partition_counter % vectorization_length;
1809  n_macro_cells_before++;
1810  }
1811  }
1812  }
1813  cell_partition_data.push_back(n_macro_cells_before);
1814  partition_l2++;
1815  }
1816  neighbor_list = neighbor_neighbor_list;
1817  neighbor_neighbor_list.resize(0);
1818  }
1819  partition_row_index[part + 1] =
1820  partition_row_index[part] + partition_l2;
1821  }
1822  }
1823  if (hp_bool == true)
1824  {
1825  partition_partition_list.swap(renumbering);
1826  }
1827  }
1828 
1829 
1830 
1831  // Function to create coloring on the second layer within each partition.
1832  // Version assumes preblocking.
1833  void
1835  const DynamicSparsityPattern & connectivity,
1836  const unsigned int partition,
1837  const std::vector<unsigned int> &cell_partition,
1838  const std::vector<unsigned int> &partition_list,
1839  const std::vector<unsigned int> &partition_size,
1840  std::vector<unsigned int> & partition_color_list)
1841  {
1842  const unsigned int n_macro_cells = *(cell_partition_data.end() - 2);
1843  std::vector<unsigned int> cell_color(n_blocks, n_macro_cells);
1844  std::vector<bool> color_finder;
1845 
1846  partition_row_index.resize(partition + 1);
1847  cell_partition_data.clear();
1848  unsigned int color_counter = 0, index_counter = 0;
1849  for (unsigned int part = 0; part < partition; part++)
1850  {
1851  partition_row_index[part] = index_counter;
1852  unsigned int max_color = 0;
1853  for (unsigned int k = partition_size[part];
1854  k < partition_size[part + 1];
1855  k++)
1856  {
1857  unsigned int cell = partition_list[k];
1858  unsigned int n_neighbors = connectivity.row_length(cell);
1859 
1860  // In the worst case, each neighbor has a different color. So we
1861  // find at least one available color between 0 and n_neighbors.
1862  color_finder.resize(n_neighbors + 1);
1863  for (unsigned int j = 0; j <= n_neighbors; ++j)
1864  color_finder[j] = true;
1866  connectivity.begin(cell),
1867  end = connectivity.end(cell);
1868  for (; neighbor != end; ++neighbor)
1869  {
1870  // Mark the color that a neighbor within the partition has
1871  // as taken
1872  if (cell_partition[neighbor->column()] == part &&
1873  cell_color[neighbor->column()] <= n_neighbors)
1874  color_finder[cell_color[neighbor->column()]] = false;
1875  }
1876  // Choose the smallest color that is not taken for the block
1877  cell_color[cell] = 0;
1878  while (color_finder[cell_color[cell]] == false)
1879  cell_color[cell]++;
1880  if (cell_color[cell] > max_color)
1881  max_color = cell_color[cell];
1882  }
1883  // Reorder within partition: First, all blocks that belong the 0 and
1884  // then so on until those with color max (Note that the smaller the
1885  // number the larger the partition)
1886  for (unsigned int color = 0; color <= max_color; color++)
1887  {
1888  cell_partition_data.push_back(color_counter);
1889  index_counter++;
1890  for (unsigned int k = partition_size[part];
1891  k < partition_size[part + 1];
1892  k++)
1893  {
1894  unsigned int cell = partition_list[k];
1895  if (cell_color[cell] == color)
1896  {
1897  partition_color_list[color_counter++] = cell;
1898  }
1899  }
1900  }
1901  }
1902  cell_partition_data.push_back(n_blocks);
1903  partition_row_index[partition] = index_counter;
1904  AssertDimension(color_counter, n_blocks);
1905  }
1906 
1907 
1908  // Function to create partitioning on the first layer.
1909  void
1911  const unsigned int cluster_size,
1912  std::vector<unsigned int> & cell_partition,
1913  std::vector<unsigned int> & partition_list,
1914  std::vector<unsigned int> & partition_size,
1915  unsigned int & partition) const
1916 
1917  {
1918  // For each block of cells, this variable saves to which partitions the
1919  // block belongs. Initialize all to n_macro_cells to mark them as not
1920  // yet assigned a partition.
1921  // std::vector<unsigned int> cell_partition (n_active_cells,
1922  // numbers::invalid_unsigned_int);
1923  // List of cells in previous partition
1924  std::vector<unsigned int> neighbor_list;
1925  // List of cells in current partition for use as neighbors in next
1926  // partition
1927  std::vector<unsigned int> neighbor_neighbor_list;
1928 
1929  // In element j of this variable, one puts the old number of the block
1930  // that should be the jth block in the new numeration.
1931  // std::vector<unsigned int> partition_list(n_active_cells,0);
1932 
1933  // This vector points to the start of each partition.
1934  // std::vector<unsigned int> partition_size(2,0);
1935 
1936  partition = 0;
1937  unsigned int counter = 0;
1938  unsigned int start_nonboundary =
1939  cell_partition_data.size() == 5 ?
1942  0;
1943 
1944  const unsigned int n_macro_cells = *(cell_partition_data.end() - 2);
1945  if (n_macro_cells == 0)
1946  return;
1947  if (scheme == color)
1948  start_nonboundary = n_macro_cells;
1949  if (scheme == partition_color ||
1950  scheme == color) // blocking_connectivity == true
1951  start_nonboundary = ((start_nonboundary + block_size - 1) / block_size);
1952  unsigned int n_blocks;
1953  if (scheme == partition_color ||
1954  scheme == color) // blocking_connectivity == true
1955  n_blocks = this->n_blocks;
1956  else
1957  n_blocks = n_active_cells;
1958 
1959  if (start_nonboundary > n_blocks)
1960  start_nonboundary = n_blocks;
1961 
1962 
1963  unsigned int start_up = 0;
1964  bool work = true;
1965  unsigned int remainder = cluster_size;
1966 
1967  // this performs a classical breath-first search in the connectivity
1968  // graph of the cells under the restriction that the size of the
1969  // partitions should be a multiple of the given block size
1970  while (work)
1971  {
1972  // put the cells with neighbors on remote MPI processes up front
1973  if (start_nonboundary > 0)
1974  {
1975  for (unsigned int cell = 0; cell < start_nonboundary; ++cell)
1976  {
1977  const unsigned int cell_nn = cell;
1978  cell_partition[cell_nn] = partition;
1979  neighbor_list.push_back(cell_nn);
1980  partition_list[counter++] = cell_nn;
1981  partition_size.back()++;
1982  }
1983  start_nonboundary = 0;
1984  remainder -= (start_nonboundary % cluster_size);
1985  if (remainder == cluster_size)
1986  remainder = 0;
1987  }
1988  else
1989  {
1990  // To start up, set the start_up cell to partition and list all
1991  // its neighbors.
1992  cell_partition[start_up] = partition;
1993  neighbor_list.push_back(start_up);
1994  partition_list[counter++] = start_up;
1995  partition_size.back()++;
1996  start_up++;
1997  remainder--;
1998  if (remainder == cluster_size)
1999  remainder = 0;
2000  }
2001  int index_before = neighbor_list.size(), index = index_before,
2002  index_stop = 0;
2003  while (remainder > 0)
2004  {
2005  if (index == index_stop)
2006  {
2007  index = neighbor_list.size();
2008  if (index == index_before)
2009  {
2010  neighbor_list.resize(0);
2011  goto not_connect;
2012  }
2013  index_stop = index_before;
2014  index_before = index;
2015  }
2016  index--;
2017  unsigned int additional = neighbor_list[index];
2019  connectivity.begin(additional),
2020  end =
2021  connectivity.end(additional);
2022  for (; neighbor != end; ++neighbor)
2023  {
2024  if (cell_partition[neighbor->column()] ==
2026  {
2027  partition_size.back()++;
2028  cell_partition[neighbor->column()] = partition;
2029  neighbor_list.push_back(neighbor->column());
2030  partition_list[counter++] = neighbor->column();
2031  remainder--;
2032  if (remainder == 0)
2033  break;
2034  }
2035  }
2036  }
2037 
2038  while (neighbor_list.size() > 0)
2039  {
2040  partition++;
2041 
2042  // counter for number of cells so far in current partition
2043  unsigned int partition_counter = 0;
2044 
2045  // Mark the start of the new partition
2046  partition_size.push_back(partition_size.back());
2047 
2048  // Loop through the list of cells in previous partition and put
2049  // all their neighbors in current partition
2050  for (const unsigned int cell : neighbor_list)
2051  {
2052  Assert(cell_partition[cell] == partition - 1,
2053  ExcInternalError());
2054  auto neighbor = connectivity.begin(cell);
2055  const auto end = connectivity.end(cell);
2056  for (; neighbor != end; ++neighbor)
2057  {
2058  if (cell_partition[neighbor->column()] ==
2060  {
2061  partition_size.back()++;
2062  cell_partition[neighbor->column()] = partition;
2063 
2064  // collect the cells of the current partition for
2065  // use as neighbors in next partition
2066  neighbor_neighbor_list.push_back(neighbor->column());
2067  partition_list[counter++] = neighbor->column();
2068  partition_counter++;
2069  }
2070  }
2071  }
2072  remainder = cluster_size - (partition_counter % cluster_size);
2073  if (remainder == cluster_size)
2074  remainder = 0;
2075  int index_stop = 0;
2076  int index_before = neighbor_neighbor_list.size(),
2077  index = index_before;
2078  while (remainder > 0)
2079  {
2080  if (index == index_stop)
2081  {
2082  index = neighbor_neighbor_list.size();
2083  if (index == index_before)
2084  {
2085  neighbor_neighbor_list.resize(0);
2086  break;
2087  }
2088  index_stop = index_before;
2089  index_before = index;
2090  }
2091  index--;
2092  unsigned int additional = neighbor_neighbor_list[index];
2094  connectivity.begin(
2095  additional),
2096  end = connectivity.end(
2097  additional);
2098  for (; neighbor != end; ++neighbor)
2099  {
2100  if (cell_partition[neighbor->column()] ==
2102  {
2103  partition_size.back()++;
2104  cell_partition[neighbor->column()] = partition;
2105  neighbor_neighbor_list.push_back(neighbor->column());
2106  partition_list[counter++] = neighbor->column();
2107  remainder--;
2108  if (remainder == 0)
2109  break;
2110  }
2111  }
2112  }
2113 
2114  neighbor_list = neighbor_neighbor_list;
2115  neighbor_neighbor_list.resize(0);
2116  }
2117  not_connect:
2118  // One has to check if the graph is not connected so we have to find
2119  // another partition.
2120  work = false;
2121  for (unsigned int j = start_up; j < n_blocks; ++j)
2122  if (cell_partition[j] == numbers::invalid_unsigned_int)
2123  {
2124  start_up = j;
2125  work = true;
2126  if (remainder == 0)
2127  remainder = cluster_size;
2128  break;
2129  }
2130  }
2131  if (remainder != 0)
2132  partition++;
2133 
2134  AssertDimension(partition_size[partition], n_blocks);
2135  }
2136 
2137 
2138  void
2139  TaskInfo::update_task_info(const unsigned int partition)
2140  {
2141  evens = (partition + 1) / 2;
2142  odds = partition / 2;
2143  n_blocked_workers = odds - (odds + evens + 1) % 2;
2145  // From here only used for partition partition option.
2146  partition_evens.resize(partition);
2147  partition_odds.resize(partition);
2148  partition_n_blocked_workers.resize(partition);
2149  partition_n_workers.resize(partition);
2150  for (unsigned int part = 0; part < partition; part++)
2151  {
2152  partition_evens[part] =
2153  (partition_row_index[part + 1] - partition_row_index[part] + 1) / 2;
2154  partition_odds[part] =
2155  (partition_row_index[part + 1] - partition_row_index[part]) / 2;
2157  partition_odds[part] -
2158  (partition_odds[part] + partition_evens[part] + 1) % 2;
2159  partition_n_workers[part] = partition_evens[part] +
2160  partition_odds[part] -
2162  }
2163  }
2164  } // namespace MatrixFreeFunctions
2165 } // namespace internal
2166 
2167 
2168 
2169 // explicit instantiations of template functions
2170 template void
2171 internal::MatrixFreeFunctions::TaskInfo::print_memory_statistics<std::ostream>(
2172  std::ostream &,
2173  const std::size_t) const;
2174 template void
2176  ConditionalOStream>(ConditionalOStream &, const std::size_t) const;
2177 
2178 
2179 DEAL_II_NAMESPACE_CLOSE
static const unsigned int invalid_unsigned_int
Definition: types.h:173
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1567
virtual void cell(const std::pair< unsigned int, unsigned int > &cell_range)=0
void parallel_for(Iterator x_begin, Iterator x_end, const Functor &functor, const unsigned int grainsize)
Definition: parallel.h:156
void loop(MFWorkerInterface &worker) const
Definition: task_info.cc:337
void add(const size_type i, const size_type j)
void guess_block_size(const unsigned int dofs_per_cell)
Definition: task_info.cc:1031
std::vector< unsigned int > partition_row_index
Definition: task_info.h:457
#define AssertIndexRange(index, range)
Definition: exceptions.h:1637
std::vector< unsigned int > partition_n_workers
Definition: task_info.h:528
virtual void vector_update_ghosts_finish()=0
Finishes the communication for the update ghost values operation.
std::vector< unsigned int > cell_partition_data
Definition: task_info.h:465
#define AssertThrow(cond, exc)
Definition: exceptions.h:1519
void make_connectivity_cells_to_blocks(const std::vector< unsigned char > &irregular_cells, const DynamicSparsityPattern &connectivity_cells, DynamicSparsityPattern &connectivity_blocks) const
Definition: task_info.cc:1469
std::vector< unsigned int > partition_n_blocked_workers
Definition: task_info.h:522
void make_partitioning_within_partitions_post_blocked(const DynamicSparsityPattern &connectivity, const std::vector< unsigned int > &cell_active_fe_index, const unsigned int partition, const unsigned int cluster_size, const bool hp_bool, const std::vector< unsigned int > &cell_partition, const std::vector< unsigned int > &partition_list, const std::vector< unsigned int > &partition_size, std::vector< unsigned int > &partition_partition_list, std::vector< unsigned char > &irregular_cells)
Definition: task_info.cc:1513
void make_coloring_within_partitions_pre_blocked(const DynamicSparsityPattern &connectivity, const unsigned int partition, const std::vector< unsigned int > &cell_partition, const std::vector< unsigned int > &partition_list, const std::vector< unsigned int > &partition_size, std::vector< unsigned int > &partition_color_list)
Definition: task_info.cc:1834
void partition(const SparsityPattern &sparsity_pattern, const unsigned int n_partitions, std::vector< unsigned int > &partition_indices, const Partitioner partitioner=Partitioner::metis)
#define Assert(cond, exc)
Definition: exceptions.h:1407
void make_thread_graph_partition_color(DynamicSparsityPattern &connectivity, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &irregular_cells, const bool hp_bool)
Definition: task_info.cc:1059
void update_task_info(const unsigned int partition)
Definition: task_info.cc:2139
std::vector< unsigned int > boundary_partition_data
Definition: task_info.h:483
void initial_setup_blocks_tasks(const std::vector< unsigned int > &boundary_cells, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &incompletely_filled_vectorization)
Definition: task_info.cc:965
virtual void zero_dst_vector_range(const unsigned int range_index)=0
virtual void boundary(const std::pair< unsigned int, unsigned int > &face_range)=0
void create_blocks_serial(const std::vector< unsigned int > &boundary_cells, const std::vector< unsigned int > &cells_close_to_boundary, const unsigned int dofs_per_cell, const std::vector< unsigned int > &cell_vectorization_categories, const bool cell_vectorization_categories_strict, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &incompletely_filled_vectorization)
Definition: task_info.cc:766
std::vector< unsigned int > invert_permutation(const std::vector< unsigned int > &permutation)
Definition: utilities.cc:790
virtual void vector_update_ghosts_start()=0
Starts the communication for the update ghost values operation.
virtual void face(const std::pair< unsigned int, unsigned int > &face_range)=0
std::vector< unsigned int > face_partition_data
Definition: task_info.h:474
void collect_boundary_cells(const unsigned int n_active_cells, const unsigned int n_active_and_ghost_cells, const unsigned int vectorization_length, std::vector< unsigned int > &boundary_cells)
Definition: task_info.cc:687
size_type row_length(const size_type row) const
void make_thread_graph_partition_partition(const std::vector< unsigned int > &cell_active_fe_index, DynamicSparsityPattern &connectivity, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &irregular_cells, const bool hp_bool)
Definition: task_info.cc:1401
void make_partitioning(const DynamicSparsityPattern &connectivity, const unsigned int cluster_size, std::vector< unsigned int > &cell_partition, std::vector< unsigned int > &partition_list, std::vector< unsigned int > &partition_size, unsigned int &partition) const
Definition: task_info.cc:1910
virtual void vector_compress_start()=0
Starts the communication for the vector compress operation.
std::vector< unsigned int > partition_odds
Definition: task_info.h:516
std::size_t memory_consumption() const
Definition: task_info.cc:670
static ::ExceptionBase & ExcNotImplemented()
static unsigned int n_threads()
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
Definition: mpi.cc:429
std::vector< unsigned int > partition_evens
Definition: task_info.h:510
void print_memory_statistics(StreamType &out, std::size_t data_length) const
Definition: task_info.cc:655
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static ::ExceptionBase & ExcInternalError()
void make_thread_graph(const std::vector< unsigned int > &cell_active_fe_index, DynamicSparsityPattern &connectivity, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &irregular_cells, const bool hp_bool)
Definition: task_info.cc:1206
virtual void vector_compress_finish()=0
Finishes the communication for the vector compress operation.