27#ifdef DEAL_II_WITH_TBB
28# include <tbb/blocked_range.h>
29# include <tbb/parallel_for.h>
31# ifndef DEAL_II_TBB_WITH_ONEAPI
32# include <tbb/task_scheduler_init.h>
62 namespace MatrixFreeFunctions
64#if defined(DEAL_II_WITH_TBB) && !defined(DEAL_II_TBB_WITH_ONEAPI)
74 ActualCellWork(MFWorkerInterface **worker_pointer,
76 const TaskInfo & task_info)
78 , worker_pointer(worker_pointer)
80 , task_info(task_info)
83 ActualCellWork(MFWorkerInterface &worker,
85 const TaskInfo & task_info)
87 , worker_pointer(nullptr)
89 , task_info(task_info)
95 MFWorkerInterface *used_worker =
96 worker !=
nullptr ? worker : *worker_pointer;
100 if (task_info.face_partition_data.empty() ==
false)
108 MFWorkerInterface * worker;
109 MFWorkerInterface **worker_pointer;
111 const TaskInfo & task_info;
114 class CellWork :
public tbb::task
117 CellWork(MFWorkerInterface &worker,
119 const TaskInfo & task_info,
120 const bool is_blocked)
123 , is_blocked(is_blocked)
131 if (is_blocked ==
true)
132 tbb::empty_task::spawn(*dummy);
136 tbb::empty_task *dummy;
140 const bool is_blocked;
145 class PartitionWork :
public tbb::task
148 PartitionWork(MFWorkerInterface &function_in,
149 const unsigned int partition_in,
150 const TaskInfo & task_info_in,
151 const bool is_blocked_in =
false)
153 , function(function_in)
155 , task_info(task_info_in)
156 , is_blocked(is_blocked_in)
162 tbb::empty_task *root =
163 new (tbb::task::allocate_root()) tbb::empty_task;
164 const unsigned int evens = task_info.partition_evens[
partition];
165 const unsigned int odds = task_info.partition_odds[
partition];
166 const unsigned int n_blocked_workers =
167 task_info.partition_n_blocked_workers[
partition];
168 const unsigned int n_workers =
169 task_info.partition_n_workers[
partition];
170 std::vector<CellWork *> worker(n_workers);
171 std::vector<CellWork *> blocked_worker(n_blocked_workers);
173 root->set_ref_count(evens + 1);
174 for (
unsigned int j = 0; j < evens; j++)
176 worker[j] =
new (root->allocate_child())
178 task_info.partition_row_index[
partition] + 2 * j,
183 worker[j]->set_ref_count(2);
184 blocked_worker[j - 1]->dummy =
185 new (worker[j]->allocate_child()) tbb::empty_task;
186 tbb::task::spawn(*blocked_worker[j - 1]);
189 worker[j]->set_ref_count(1);
192 blocked_worker[j] =
new (worker[j]->allocate_child())
194 task_info.partition_row_index[
partition] + 2 * j +
203 worker[evens] =
new (worker[j]->allocate_child())
205 task_info.partition_row_index[
partition] +
209 tbb::task::spawn(*worker[evens]);
213 tbb::empty_task *child =
214 new (worker[j]->allocate_child()) tbb::empty_task();
215 tbb::task::spawn(*child);
220 root->wait_for_all();
221 root->destroy(*root);
222 if (is_blocked ==
true)
223 tbb::empty_task::spawn(*dummy);
227 tbb::empty_task *dummy;
230 MFWorkerInterface &function;
232 const TaskInfo & task_info;
233 const bool is_blocked;
245 CellWork(MFWorkerInterface &worker_in,
246 const TaskInfo & task_info_in,
247 const unsigned int partition_in)
249 , task_info(task_info_in)
254 operator()(
const tbb::blocked_range<unsigned int> &r)
const
256 const unsigned int start_index =
257 task_info.cell_partition_data[
partition] +
258 task_info.block_size * r.begin();
259 const unsigned int end_index =
260 std::min(start_index + task_info.block_size * (r.end() - r.begin()),
261 task_info.cell_partition_data[
partition + 1]);
262 worker.cell(std::make_pair(start_index, end_index));
264 if (task_info.face_partition_data.empty() ==
false)
271 MFWorkerInterface &worker;
272 const TaskInfo & task_info;
278 class PartitionWork :
public tbb::task
281 PartitionWork(MFWorkerInterface &worker_in,
282 const unsigned int partition_in,
283 const TaskInfo & task_info_in,
284 const bool is_blocked_in)
288 , task_info(task_info_in)
289 , is_blocked(is_blocked_in)
295 const unsigned int n_chunks =
296 (task_info.cell_partition_data[
partition + 1] -
297 task_info.cell_partition_data[
partition] + task_info.block_size -
299 task_info.block_size;
300 parallel_for(tbb::blocked_range<unsigned int>(0, n_chunks, 1),
302 if (is_blocked ==
true)
303 tbb::empty_task::spawn(*dummy);
307 tbb::empty_task *dummy;
310 MFWorkerInterface &worker;
312 const TaskInfo & task_info;
313 const bool is_blocked;
320 class MPICommunication :
public tbb::task
323 MPICommunication(MFWorkerInterface &worker_in,
const bool do_compress)
325 , do_compress(do_compress)
331 if (do_compress ==
false)
332 worker.vector_update_ghosts_finish();
334 worker.vector_compress_start();
339 MFWorkerInterface &worker;
340 const bool do_compress;
362#if defined(DEAL_II_WITH_TBB) && !defined(DEAL_II_TBB_WITH_ONEAPI)
369 tbb::empty_task *root =
370 new (tbb::task::allocate_root()) tbb::empty_task;
371 root->set_ref_count(
evens + 1);
372 std::vector<partition::PartitionWork *> worker(
n_workers);
373 std::vector<partition::PartitionWork *> blocked_worker(
375 MPICommunication *worker_compr =
376 new (root->allocate_child()) MPICommunication(funct,
true);
377 worker_compr->set_ref_count(1);
378 for (
unsigned int j = 0; j <
evens; j++)
382 worker[j] =
new (root->allocate_child())
383 partition::PartitionWork(funct, 2 * j, *
this,
false);
384 worker[j]->set_ref_count(2);
385 blocked_worker[j - 1]->dummy =
386 new (worker[j]->allocate_child()) tbb::empty_task;
387 tbb::task::spawn(*blocked_worker[j - 1]);
391 worker[j] =
new (worker_compr->allocate_child())
392 partition::PartitionWork(funct, 2 * j, *
this,
false);
393 worker[j]->set_ref_count(2);
394 MPICommunication *worker_dist =
395 new (worker[j]->allocate_child())
396 MPICommunication(funct,
false);
397 tbb::task::spawn(*worker_dist);
401 blocked_worker[j] =
new (worker[j]->allocate_child())
402 partition::PartitionWork(funct, 2 * j + 1, *
this,
true);
408 worker[
evens] =
new (worker[j]->allocate_child())
409 partition::PartitionWork(funct,
413 tbb::task::spawn(*worker[
evens]);
417 tbb::empty_task *child =
418 new (worker[j]->allocate_child()) tbb::empty_task();
419 tbb::task::spawn(*child);
424 root->wait_for_all();
425 root->destroy(*root);
441 tbb::empty_task *root =
442 new (tbb::task::allocate_root()) tbb::empty_task;
443 root->set_ref_count(
evens + 1);
448 std::vector<color::PartitionWork *> worker(
n_workers);
449 std::vector<color::PartitionWork *> blocked_worker(
451 unsigned int worker_index = 0, slice_index = 0;
452 int spawn_index_child = -2;
453 MPICommunication *worker_compr =
454 new (root->allocate_child()) MPICommunication(funct,
true);
455 worker_compr->set_ref_count(1);
456 for (
unsigned int part = 0;
461 worker[worker_index] =
462 new (worker_compr->allocate_child())
463 color::PartitionWork(funct,
468 worker[worker_index] =
new (root->allocate_child())
469 color::PartitionWork(funct,
477 worker[worker_index]->set_ref_count(1);
479 worker[worker_index] =
480 new (worker[worker_index - 1]->allocate_child())
481 color::PartitionWork(funct,
486 worker[worker_index]->set_ref_count(2);
489 blocked_worker[(part - 1) / 2]->dummy =
490 new (worker[worker_index]->allocate_child())
493 if (spawn_index_child == -1)
494 tbb::task::spawn(*blocked_worker[(part - 1) / 2]);
497 Assert(spawn_index_child >= 0,
499 tbb::task::spawn(*worker[spawn_index_child]);
501 spawn_index_child = -2;
505 MPICommunication *worker_dist =
506 new (worker[worker_index]->allocate_child())
507 MPICommunication(funct,
false);
508 tbb::task::spawn(*worker_dist);
516 blocked_worker[part / 2] =
517 new (worker[worker_index - 1]->allocate_child())
518 color::PartitionWork(funct,
525 blocked_worker[part / 2]->set_ref_count(1);
526 worker[worker_index] =
new (
527 blocked_worker[part / 2]->allocate_child())
528 color::PartitionWork(funct,
536 spawn_index_child = -1;
545 worker[worker_index]->set_ref_count(1);
548 worker[worker_index] =
549 new (worker[worker_index - 1]->allocate_child())
550 color::PartitionWork(funct,
555 spawn_index_child = worker_index;
560 tbb::empty_task *
final =
561 new (worker[worker_index - 1]->allocate_child())
563 tbb::task::spawn(*
final);
564 spawn_index_child = worker_index - 1;
570 tbb::task::spawn(*worker[spawn_index_child]);
572 root->wait_for_all();
573 root->destroy(*root);
586 tbb::empty_task *root =
587 new (tbb::task::allocate_root()) tbb::empty_task;
588 root->set_ref_count(2);
589 color::PartitionWork *worker =
590 new (root->allocate_child())
591 color::PartitionWork(funct,
color, *
this,
false);
592 tbb::empty_task::spawn(*worker);
593 root->wait_for_all();
594 root->destroy(*root);
686 template <
typename StreamType>
689 const std::size_t data_length)
const
696 out << memory_c.
min <<
"/" << memory_c.
avg <<
"/" << memory_c.
max;
697 out <<
" MB" << std::endl;
721 std::vector<unsigned int> &boundary_cells)
725 unsigned int fillup_needed =
733 std::vector<unsigned int> new_boundary_cells;
734 new_boundary_cells.reserve(boundary_cells.size());
736 unsigned int next_free_slot = 0, bound_index = 0;
737 while (fillup_needed > 0 && bound_index < boundary_cells.size())
739 if (next_free_slot < boundary_cells[bound_index])
743 if (next_free_slot + fillup_needed <=
744 boundary_cells[bound_index])
746 for (
unsigned int j =
747 boundary_cells[bound_index] - fillup_needed;
748 j < boundary_cells[bound_index];
750 new_boundary_cells.push_back(j);
757 for (
unsigned int j = next_free_slot;
758 j < boundary_cells[bound_index];
760 new_boundary_cells.push_back(j);
762 boundary_cells[bound_index] - next_free_slot;
765 new_boundary_cells.push_back(boundary_cells[bound_index]);
766 next_free_slot = boundary_cells[bound_index] + 1;
769 while (fillup_needed > 0 &&
770 (new_boundary_cells.size() == 0 ||
772 new_boundary_cells.push_back(new_boundary_cells.back() + 1);
773 while (bound_index < boundary_cells.size())
774 new_boundary_cells.push_back(boundary_cells[bound_index++]);
776 boundary_cells.swap(new_boundary_cells);
780 std::sort(boundary_cells.begin(), boundary_cells.end());
793 const std::vector<unsigned int> &cells_with_comm,
794 const unsigned int dofs_per_cell,
795 const bool categories_are_hp,
796 const std::vector<unsigned int> &cell_vectorization_categories,
797 const bool cell_vectorization_categories_strict,
798 const std::vector<unsigned int> &parent_relation,
799 std::vector<unsigned int> & renumbering,
800 std::vector<unsigned char> & incompletely_filled_vectorization)
826 unsigned int vectorization_length_bits = 0;
828 while (my_length >>= 1)
829 ++vectorization_length_bits;
830 const unsigned int n_lanes = 1 << vectorization_length_bits;
835 unsigned int n_categories = 1;
837 if (cell_vectorization_categories.empty() ==
false)
842 std::set<unsigned int> used_categories;
844 used_categories.insert(cell_vectorization_categories[i]);
845 std::vector<unsigned int> used_categories_vector(
846 used_categories.size());
848 for (
const auto &it : used_categories)
849 used_categories_vector[n_categories++] = it;
852 const unsigned int index =
854 used_categories_vector.end(),
855 cell_vectorization_categories[i]) -
856 used_categories_vector.begin();
858 tight_category_map[i] = index;
865 std::vector<std::vector<unsigned int>> renumbering_category(n_categories);
867 renumbering_category[tight_category_map[i]].push_back(i);
869 if (cell_vectorization_categories_strict ==
false && n_categories > 1)
870 for (
unsigned int j = n_categories - 1; j > 0; --j)
872 unsigned int lower_index = j - 1;
873 while (renumbering_category[j].size() % n_lanes)
875 while (renumbering_category[j].size() % n_lanes &&
876 !renumbering_category[lower_index].empty())
878 renumbering_category[j].push_back(
879 renumbering_category[lower_index].back());
880 renumbering_category[lower_index].pop_back();
882 if (lower_index == 0)
893 std::vector<unsigned int> temporary_numbering;
895 (n_lanes - 1) * n_categories);
896 const unsigned int n_cells_per_parent =
897 std::count(parent_relation.begin(), parent_relation.end(), 0);
898 std::vector<unsigned int> category_size;
899 for (
unsigned int j = 0; j < n_categories; ++j)
901 std::vector<std::pair<unsigned int, unsigned int>> grouped_cells;
902 std::vector<unsigned int> other_cells;
903 for (
const unsigned int cell : renumbering_category[j])
904 if (parent_relation.empty() ||
906 other_cells.push_back(cell);
908 grouped_cells.emplace_back(parent_relation[cell], cell);
911 std::sort(grouped_cells.begin(), grouped_cells.end());
912 std::vector<unsigned int> n_cells_per_group;
913 unsigned int length = 0;
914 for (
unsigned int i = 0; i < grouped_cells.size(); ++i, ++length)
915 if (i > 0 && grouped_cells[i].
first != grouped_cells[i - 1].
first)
917 n_cells_per_group.push_back(length);
921 n_cells_per_group.push_back(length);
926 auto group_it = grouped_cells.begin();
927 for (
unsigned int length : n_cells_per_group)
928 if (length < n_cells_per_parent)
929 for (
unsigned int j = 0; j < length; ++j)
930 other_cells.push_back((group_it++)->second);
936 for (
unsigned int j = 0; j < length; ++j)
937 temporary_numbering.push_back((group_it++)->second);
941 std::sort(other_cells.begin(), other_cells.end());
942 temporary_numbering.insert(temporary_numbering.end(),
946 while (temporary_numbering.size() % n_lanes != 0)
949 category_size.push_back(temporary_numbering.size());
953 std::vector<bool> batch_with_comm(temporary_numbering.size() / n_lanes,
955 std::vector<unsigned int> temporary_numbering_inverse(
n_active_cells);
956 for (
unsigned int i = 0; i < temporary_numbering.size(); ++i)
958 temporary_numbering_inverse[temporary_numbering[i]] = i;
959 for (
const unsigned int cell : cells_with_comm)
960 batch_with_comm[temporary_numbering_inverse[cell] / n_lanes] =
true;
966 std::vector<std::array<unsigned int, 3>> batch_order;
967 std::vector<std::array<unsigned int, 3>> batch_order_comm;
968 for (
unsigned int i = 0; i < temporary_numbering.size(); i += n_lanes)
970 unsigned int max_index = 0;
971 for (
unsigned int j = 0; j < n_lanes; ++j)
973 max_index =
std::max(temporary_numbering[i + j], max_index);
974 const unsigned int category_hp =
976 std::upper_bound(category_size.begin(), category_size.end(), i) -
977 category_size.begin() :
979 const std::array<unsigned int, 3> next{{category_hp, max_index, i}};
980 if (batch_with_comm[i / n_lanes])
981 batch_order_comm.emplace_back(next);
983 batch_order.emplace_back(next);
986 std::sort(batch_order.begin(), batch_order.end());
987 std::sort(batch_order_comm.begin(), batch_order_comm.end());
994 std::vector<unsigned int> blocks;
997 if (batch_order.empty())
998 std::swap(batch_order_comm, batch_order);
1001 blocks = {0,
static_cast<unsigned int>(batch_order.size())};
1006 const unsigned int comm_begin = batch_order.size() / 2;
1007 batch_order.insert(batch_order.begin() + comm_begin,
1008 batch_order_comm.begin(),
1009 batch_order_comm.end());
1010 const unsigned int comm_end = comm_begin + batch_order_comm.size();
1011 const unsigned int end = batch_order.size();
1012 blocks = {0, comm_begin, comm_end,
end};
1016 const unsigned int n_cell_batches = batch_order.size();
1017 const unsigned int n_ghost_batches =
1019 incompletely_filled_vectorization.resize(n_cell_batches +
1025 renumbering.clear();
1029 unsigned int counter = 0;
1030 for (
unsigned int block = 0; block < blocks.size() - 1; ++block)
1032 const unsigned int grain_size =
1033 std::max((2048U / dofs_per_cell) / 8 * 4, 2U);
1034 for (
unsigned int k = blocks[block]; k < blocks[block + 1];
1037 std::min(k + grain_size, blocks[block + 1]));
1041 for (
unsigned int k = blocks[block]; k < blocks[block + 1]; ++k)
1043 const unsigned int pos = batch_order[k][2];
1045 for (; j < n_lanes && temporary_numbering[pos + j] !=
1048 renumbering[counter++] = temporary_numbering[pos + j];
1050 incompletely_filled_vectorization[k] = j;
1060 if (!cell_vectorization_categories.empty())
1063 renumbering[cell] = cell;
1066 incompletely_filled_vectorization.back() =
n_ghost_cells % n_lanes;
1071 std::vector<unsigned int> renumber_cpy(renumbering);
1072 std::sort(renumber_cpy.begin(), renumber_cpy.end());
1073 for (
unsigned int i = 0; i < renumber_cpy.size(); ++i)
1082 const std::vector<unsigned int> &boundary_cells,
1083 std::vector<unsigned int> & renumbering,
1084 std::vector<unsigned char> & incompletely_filled_vectorization)
1086 const unsigned int n_cell_batches =
1088 const unsigned int n_ghost_slots =
1090 incompletely_filled_vectorization.resize(n_cell_batches + n_ghost_slots);
1092 incompletely_filled_vectorization[n_cell_batches - 1] =
1096 incompletely_filled_vectorization[n_cell_batches + n_ghost_slots - 1] =
1100 std::vector<unsigned int> reverse_numbering(
1102 for (
unsigned int j = 0; j < boundary_cells.size(); ++j)
1103 reverse_numbering[boundary_cells[j]] = j;
1104 unsigned int counter = boundary_cells.size();
1107 reverse_numbering[j] = counter++;
1114 renumbering.push_back(j);
1122 const unsigned int n_macro_boundary_cells =
1126 (n_cell_batches - n_macro_boundary_cells) / 2);
1128 n_macro_boundary_cells);
1166 1 <<
static_cast<unsigned int>(std::log2(
block_size + 1));
1177 std::vector<unsigned int> & renumbering,
1178 std::vector<unsigned char> &irregular_cells,
1182 if (n_cell_batches == 0)
1187 unsigned int partition = 0, counter = 0;
1201 std::vector<unsigned int> cell_partition(
n_blocks,
1206 std::vector<unsigned int> partition_list(
n_blocks, 0);
1207 std::vector<unsigned int> partition_color_list(
n_blocks, 0);
1210 std::vector<unsigned int> partition_size(2, 0);
1216 unsigned int cluster_size = 1;
1232 partition_color_list);
1234 partition_list = renumbering;
1239 std::vector<unsigned int> sorted_pc_list(partition_color_list);
1240 std::sort(sorted_pc_list.begin(), sorted_pc_list.end());
1241 for (
unsigned int i = 0; i < sorted_pc_list.size(); ++i)
1248 std::vector<unsigned int> block_start(n_cell_batches + 1);
1249 std::vector<unsigned char> irregular(n_cell_batches);
1251 unsigned int mcell_start = 0;
1253 for (
unsigned int block = 0; block <
n_blocks; block++)
1255 block_start[block + 1] = block_start[block];
1256 for (
unsigned int mcell = mcell_start;
1260 unsigned int n_comp = (irregular_cells[mcell] > 0) ?
1261 irregular_cells[mcell] :
1263 block_start[block + 1] += n_comp;
1269 unsigned int counter_macro = 0;
1270 unsigned int block_size_last =
1272 if (block_size_last == 0)
1275 unsigned int tick = 0;
1276 for (
unsigned int block = 0; block <
n_blocks; block++)
1278 unsigned int present_block = partition_color_list[block];
1279 for (
unsigned int cell = block_start[present_block];
1280 cell < block_start[present_block + 1];
1282 renumbering[counter++] = partition_list[cell];
1283 unsigned int this_block_size =
1291 for (
unsigned int j = 0; j < this_block_size; j++)
1292 irregular[counter_macro++] =
1293 irregular_cells[present_block *
block_size + j];
1298 irregular_cells.swap(irregular);
1305 std::vector<unsigned int> sorted_renumbering(renumbering);
1306 std::sort(sorted_renumbering.begin(), sorted_renumbering.end());
1307 for (
unsigned int i = 0; i < sorted_renumbering.size(); ++i)
1323 const std::vector<unsigned int> &cell_active_fe_index,
1325 std::vector<unsigned int> & renumbering,
1326 std::vector<unsigned char> & irregular_cells,
1330 if (n_cell_batches == 0)
1340 connectivity_blocks);
1352 std::vector<unsigned int> cell_partition(
n_blocks,
1358 std::vector<unsigned int> partition_list(
n_blocks, 0);
1359 std::vector<unsigned int> partition_2layers_list(
n_blocks, 0);
1362 std::vector<unsigned int> partition_size(2, 0);
1370 unsigned int cluster_size = 1;
1397 cell_active_fe_index,
1404 partition_2layers_list,
1414 partition_2layers_list);
1420 std::vector<unsigned int> sorted_pc_list(partition_2layers_list);
1421 std::sort(sorted_pc_list.begin(), sorted_pc_list.end());
1422 for (
unsigned int i = 0; i < sorted_pc_list.size(); ++i)
1429 renumbering_in.swap(renumbering);
1435 for (
unsigned int j = 0; j < renumbering.size(); j++)
1436 renumbering[j] = renumbering_in[partition_2layers_list[j]];
1445 std::vector<unsigned int> block_start(n_cell_batches + 1);
1446 std::vector<unsigned char> irregular(n_cell_batches);
1448 unsigned int counter = 0;
1449 unsigned int mcell_start = 0;
1451 for (
unsigned int block = 0; block <
n_blocks; block++)
1453 block_start[block + 1] = block_start[block];
1454 for (
unsigned int mcell = mcell_start;
1458 unsigned int n_comp = (irregular_cells[mcell] > 0) ?
1459 irregular_cells[mcell] :
1461 block_start[block + 1] += n_comp;
1467 unsigned int counter_macro = 0;
1468 unsigned int block_size_last =
1470 if (block_size_last == 0)
1473 unsigned int tick = 0;
1474 for (
unsigned int block = 0; block <
n_blocks; block++)
1476 unsigned int present_block = partition_2layers_list[block];
1477 for (
unsigned int cell = block_start[present_block];
1478 cell < block_start[present_block + 1];
1480 renumbering[counter++] = renumbering_in[cell];
1481 unsigned int this_block_size =
1489 for (
unsigned int j = 0; j < this_block_size; j++)
1490 irregular[counter_macro++] =
1491 irregular_cells[present_block *
block_size + j];
1496 irregular_cells.swap(irregular);
1502 std::vector<unsigned int> sorted_renumbering(renumbering);
1503 std::sort(sorted_renumbering.begin(), sorted_renumbering.end());
1504 for (
unsigned int i = 0; i < sorted_renumbering.size(); ++i)
1518 const std::vector<unsigned int> &cell_active_fe_index,
1520 std::vector<unsigned int> & renumbering,
1521 std::vector<unsigned char> & irregular_cells,
1525 if (n_cell_batches == 0)
1542 std::vector<unsigned int> partition_partition_list(
n_active_cells, 0);
1545 std::vector<unsigned int> partition_size(2, 0);
1561 cell_active_fe_index,
1568 partition_partition_list,
1571 partition_list.swap(renumbering);
1573 for (
unsigned int j = 0; j < renumbering.size(); j++)
1574 renumbering[j] = partition_list[partition_partition_list[j]];
1586 const std::vector<unsigned char> &irregular_cells,
1590 std::vector<std::vector<unsigned int>> cell_blocks(
n_blocks);
1592 unsigned int cell = 0;
1593 for (
unsigned int i = 0, mcell = 0; i <
n_blocks; ++i)
1595 for (
unsigned int c = 0;
1599 unsigned int ncomp = (irregular_cells[mcell] > 0) ?
1600 irregular_cells[mcell] :
1602 for (
unsigned int c = 0; c < ncomp; ++c, ++cell)
1604 cell_blocks[i].push_back(cell);
1605 touched_cells[cell] = i;
1610 for (
unsigned int i = 0; i < cell_blocks.size(); ++i)
1611 for (
unsigned int col = 0; col < cell_blocks[i].size(); ++col)
1614 connectivity_cells.
begin(cell_blocks[i][col]);
1615 it != connectivity_cells.
end(cell_blocks[i][col]);
1618 if (touched_cells[it->column()] != i)
1619 connectivity_blocks.
add(i, touched_cells[it->column()]);
1631 const std::vector<unsigned int> &cell_active_fe_index,
1633 const unsigned int cluster_size,
1635 const std::vector<unsigned int> &cell_partition,
1636 const std::vector<unsigned int> &partition_list,
1637 const std::vector<unsigned int> &partition_size,
1638 std::vector<unsigned int> & partition_partition_list,
1639 std::vector<unsigned char> & irregular_cells)
1642 const unsigned int n_ghost_slots =
1646 std::vector<unsigned int> neighbor_list;
1649 std::vector<unsigned int> neighbor_neighbor_list;
1653 irregular_cells.back() = 0;
1656 unsigned int max_fe_index = 0;
1657 for (
const unsigned int fe_index : cell_active_fe_index)
1658 max_fe_index =
std::max(fe_index, max_fe_index);
1664 unsigned int n_cell_batches_before = 0;
1670 std::vector<unsigned int> cell_partition_l2(
1676 unsigned int counter = 0;
1677 unsigned int missing_macros;
1678 for (
unsigned int part = 0; part <
partition; ++part)
1680 neighbor_neighbor_list.resize(0);
1681 neighbor_list.resize(0);
1683 unsigned int partition_l2 = 0;
1684 unsigned int start_up = partition_size[part];
1685 unsigned int partition_counter = 0;
1688 if (neighbor_list.size() == 0)
1691 partition_counter = 0;
1692 for (
unsigned int j = start_up;
1693 j < partition_size[part + 1];
1695 if (cell_partition[partition_list[j]] == part &&
1696 cell_partition_l2[partition_list[j]] ==
1701 partition_counter = 1;
1705 cell_partition_l2[partition_list[start_up]] =
1707 neighbor_neighbor_list.push_back(
1708 partition_list[start_up]);
1709 partition_partition_list[counter++] =
1710 partition_list[start_up];
1717 partition_counter = 0;
1718 for (
const unsigned int neighbor : neighbor_list)
1720 Assert(cell_partition[neighbor] == part,
1722 Assert(cell_partition_l2[neighbor] == partition_l2 - 1,
1724 auto neighbor_it = connectivity.
begin(neighbor);
1725 const auto end_it = connectivity.
end(neighbor);
1726 for (; neighbor_it != end_it; ++neighbor_it)
1728 if (cell_partition[neighbor_it->column()] == part &&
1729 cell_partition_l2[neighbor_it->column()] ==
1732 cell_partition_l2[neighbor_it->column()] =
1734 neighbor_neighbor_list.push_back(
1735 neighbor_it->column());
1736 partition_partition_list[counter++] =
1737 neighbor_it->column();
1738 partition_counter++;
1743 if (partition_counter > 0)
1745 int index_before = neighbor_neighbor_list.size(),
1746 index = index_before;
1751 std::vector<unsigned int> remaining_per_cell_batch(
1753 std::vector<std::vector<unsigned int>>
1754 renumbering_fe_index;
1757 if (hp_bool ==
true)
1759 renumbering_fe_index.resize(max_fe_index + 1);
1760 for (cell = counter - partition_counter;
1764 renumbering_fe_index
1765 [cell_active_fe_index.empty() ?
1767 cell_active_fe_index
1768 [partition_partition_list[cell]]]
1769 .push_back(partition_partition_list[cell]);
1772 for (
unsigned int j = 0; j < max_fe_index + 1; j++)
1774 remaining_per_cell_batch[j] =
1775 renumbering_fe_index[j].size() %
1777 if (remaining_per_cell_batch[j] != 0)
1780 ((renumbering_fe_index[j].size() +
1787 remaining_per_cell_batch.resize(1);
1788 remaining_per_cell_batch[0] =
1792 if (remaining_per_cell_batch[0] != 0)
1799 cluster_size - (missing_macros % cluster_size);
1802 while (missing_macros > 0 || filled ==
false)
1806 index = neighbor_neighbor_list.size();
1807 if (index == index_before)
1809 if (missing_macros != 0)
1811 neighbor_neighbor_list.resize(0);
1816 index_before = index;
1819 unsigned int additional =
1820 neighbor_neighbor_list[index];
1831 for (; neighbor !=
end; ++neighbor)
1833 if (cell_partition[neighbor->
column()] == part &&
1834 cell_partition_l2[neighbor->
column()] ==
1837 unsigned int this_index = 0;
1838 if (hp_bool ==
true)
1840 cell_active_fe_index.empty() ?
1842 cell_active_fe_index[neighbor
1849 if (missing_macros > 0 ||
1850 remaining_per_cell_batch[this_index] > 0)
1852 cell_partition_l2[neighbor->
column()] =
1854 neighbor_neighbor_list.push_back(
1856 if (hp_bool ==
true)
1857 renumbering_fe_index[this_index]
1858 .push_back(neighbor->
column());
1859 partition_partition_list[counter] =
1862 partition_counter++;
1863 if (remaining_per_cell_batch
1864 [this_index] == 0 &&
1867 remaining_per_cell_batch[this_index]++;
1868 if (remaining_per_cell_batch
1872 remaining_per_cell_batch[this_index] =
1875 if (missing_macros == 0)
1878 for (
unsigned int fe_ind = 0;
1879 fe_ind < max_fe_index + 1;
1881 if (remaining_per_cell_batch
1891 if (hp_bool ==
true)
1896 cell = counter - partition_counter;
1897 for (
unsigned int j = 0; j < max_fe_index + 1; j++)
1899 for (
const unsigned int jj :
1900 renumbering_fe_index[j])
1901 renumbering[cell++] = jj;
1902 if (renumbering_fe_index[j].size() %
1905 irregular_cells[renumbering_fe_index[j].size() /
1907 n_cell_batches_before] =
1908 renumbering_fe_index[j].size() %
1910 n_cell_batches_before +=
1911 (renumbering_fe_index[j].size() +
1914 renumbering_fe_index[j].resize(0);
1919 n_cell_batches_before +=
1923 irregular_cells[n_cell_batches_before] =
1925 n_cell_batches_before++;
1932 neighbor_list = neighbor_neighbor_list;
1933 neighbor_neighbor_list.resize(0);
1939 if (hp_bool ==
true)
1941 partition_partition_list.swap(renumbering);
1953 const std::vector<unsigned int> &cell_partition,
1954 const std::vector<unsigned int> &partition_list,
1955 const std::vector<unsigned int> &partition_size,
1956 std::vector<unsigned int> & partition_color_list)
1959 std::vector<unsigned int> cell_color(
n_blocks, n_cell_batches);
1960 std::vector<bool> color_finder;
1964 unsigned int color_counter = 0, index_counter = 0;
1965 for (
unsigned int part = 0; part <
partition; part++)
1968 unsigned int max_color = 0;
1969 for (
unsigned int k = partition_size[part];
1970 k < partition_size[part + 1];
1973 unsigned int cell = partition_list[k];
1974 unsigned int n_neighbors = connectivity.
row_length(cell);
1978 color_finder.resize(n_neighbors + 1);
1979 for (
unsigned int j = 0; j <= n_neighbors; ++j)
1980 color_finder[j] =
true;
1982 connectivity.
begin(cell),
1983 end = connectivity.
end(cell);
1984 for (; neighbor !=
end; ++neighbor)
1988 if (cell_partition[neighbor->
column()] == part &&
1989 cell_color[neighbor->
column()] <= n_neighbors)
1990 color_finder[cell_color[neighbor->
column()]] =
false;
1993 cell_color[cell] = 0;
1994 while (color_finder[cell_color[cell]] ==
false)
1996 if (cell_color[cell] > max_color)
1997 max_color = cell_color[cell];
2006 for (
unsigned int k = partition_size[part];
2007 k < partition_size[part + 1];
2010 unsigned int cell = partition_list[k];
2011 if (cell_color[cell] ==
color)
2013 partition_color_list[color_counter++] = cell;
2027 const unsigned int cluster_size,
2028 std::vector<unsigned int> & cell_partition,
2029 std::vector<unsigned int> & partition_list,
2030 std::vector<unsigned int> & partition_size,
2040 std::vector<unsigned int> neighbor_list;
2043 std::vector<unsigned int> neighbor_neighbor_list;
2053 unsigned int counter = 0;
2054 unsigned int start_nonboundary =
2061 if (n_cell_batches == 0)
2064 start_nonboundary = n_cell_batches;
2079 unsigned int start_up = 0;
2081 unsigned int remainder = cluster_size;
2089 if (start_nonboundary > 0)
2091 for (
unsigned int cell = 0; cell < start_nonboundary; ++cell)
2093 const unsigned int cell_nn = cell;
2095 neighbor_list.push_back(cell_nn);
2096 partition_list[counter++] = cell_nn;
2097 partition_size.back()++;
2099 start_nonboundary = 0;
2100 remainder -= (start_nonboundary % cluster_size);
2101 if (remainder == cluster_size)
2109 neighbor_list.push_back(start_up);
2110 partition_list[counter++] = start_up;
2111 partition_size.back()++;
2114 if (remainder == cluster_size)
2117 int index_before = neighbor_list.size(), index = index_before,
2119 while (remainder > 0)
2121 if (index == index_stop)
2123 index = neighbor_list.size();
2124 if (index == index_before)
2126 neighbor_list.resize(0);
2129 index_stop = index_before;
2130 index_before = index;
2133 unsigned int additional = neighbor_list[index];
2135 connectivity.
begin(additional),
2137 connectivity.
end(additional);
2138 for (; neighbor !=
end; ++neighbor)
2140 if (cell_partition[neighbor->
column()] ==
2143 partition_size.back()++;
2145 neighbor_list.push_back(neighbor->
column());
2146 partition_list[counter++] = neighbor->
column();
2154 while (neighbor_list.size() > 0)
2159 unsigned int partition_counter = 0;
2162 partition_size.push_back(partition_size.back());
2166 for (
const unsigned int cell : neighbor_list)
2170 auto neighbor = connectivity.
begin(cell);
2171 const auto end = connectivity.
end(cell);
2172 for (; neighbor !=
end; ++neighbor)
2174 if (cell_partition[neighbor->column()] ==
2177 partition_size.back()++;
2178 cell_partition[neighbor->column()] =
partition;
2182 neighbor_neighbor_list.push_back(neighbor->column());
2183 partition_list[counter++] = neighbor->column();
2184 partition_counter++;
2188 remainder = cluster_size - (partition_counter % cluster_size);
2189 if (remainder == cluster_size)
2192 int index_before = neighbor_neighbor_list.size(),
2193 index = index_before;
2194 while (remainder > 0)
2196 if (index == index_stop)
2198 index = neighbor_neighbor_list.size();
2199 if (index == index_before)
2201 neighbor_neighbor_list.resize(0);
2204 index_stop = index_before;
2205 index_before = index;
2208 unsigned int additional = neighbor_neighbor_list[index];
2214 for (; neighbor !=
end; ++neighbor)
2216 if (cell_partition[neighbor->
column()] ==
2219 partition_size.back()++;
2221 neighbor_neighbor_list.push_back(neighbor->
column());
2222 partition_list[counter++] = neighbor->
column();
2230 neighbor_list = neighbor_neighbor_list;
2231 neighbor_neighbor_list.resize(0);
2237 for (
unsigned int j = start_up; j <
n_blocks; ++j)
2243 remainder = cluster_size;
2266 for (
unsigned int part = 0; part <
partition; part++)
2287internal::MatrixFreeFunctions::TaskInfo::print_memory_statistics<std::ostream>(
2289 const std::size_t)
const;
static unsigned int n_threads()
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
#define AssertThrow(cond, exc)
size_type row_length(const size_type row) const
void add(const size_type i, const size_type j)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
void swap(MemorySpaceData< Number, MemorySpace > &, MemorySpaceData< Number, MemorySpace > &)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
VectorType::value_type * end(VectorType &V)
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
Iterator lower_bound(Iterator first, Iterator last, const T &val)
std::vector< Integer > invert_permutation(const std::vector< Integer > &permutation)
unsigned int minimum_parallel_grain_size
static const unsigned int invalid_unsigned_int
void parallel_for(Iterator x_begin, Iterator x_end, const Functor &functor, const unsigned int grainsize)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
virtual void face(const unsigned int range_index)=0
virtual void zero_dst_vector_range(const unsigned int range_index)=0
virtual void cell(const std::pair< unsigned int, unsigned int > &cell_range)=0
virtual void boundary(const unsigned int range_index)=0
virtual void vector_compress_start()=0
Starts the communication for the vector compress operation.
virtual void cell_loop_post_range(const unsigned int range_index)=0
virtual void vector_update_ghosts_start()=0
Starts the communication for the update ghost values operation.
virtual void cell_loop_pre_range(const unsigned int range_index)=0
virtual void vector_update_ghosts_finish()=0
Finishes the communication for the update ghost values operation.
virtual void vector_compress_finish()=0
Finishes the communication for the vector compress operation.
unsigned int n_ghost_cells
std::size_t memory_consumption() const
std::vector< unsigned int > boundary_partition_data
void loop(MFWorkerInterface &worker) const
void make_thread_graph_partition_color(DynamicSparsityPattern &connectivity, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &irregular_cells, const bool hp_bool)
std::vector< unsigned int > partition_n_workers
unsigned int vectorization_length
void create_blocks_serial(const std::vector< unsigned int > &cells_with_comm, const unsigned int dofs_per_cell, const bool categories_are_hp, const std::vector< unsigned int > &cell_vectorization_categories, const bool cell_vectorization_categories_strict, const std::vector< unsigned int > &parent_relation, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &incompletely_filled_vectorization)
void make_connectivity_cells_to_blocks(const std::vector< unsigned char > &irregular_cells, const DynamicSparsityPattern &connectivity_cells, DynamicSparsityPattern &connectivity_blocks) const
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)
unsigned int n_active_cells
void initial_setup_blocks_tasks(const std::vector< unsigned int > &boundary_cells, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &incompletely_filled_vectorization)
void print_memory_statistics(StreamType &out, std::size_t data_length) const
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)
std::vector< unsigned int > partition_row_index
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)
std::vector< unsigned int > partition_evens
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)
unsigned int n_blocked_workers
std::vector< unsigned int > partition_n_blocked_workers
std::vector< unsigned int > cell_partition_data
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
void update_task_info(const unsigned int partition)
void make_boundary_cells_divisible(std::vector< unsigned int > &boundary_cells)
TasksParallelScheme scheme
void guess_block_size(const unsigned int dofs_per_cell)
std::vector< unsigned int > partition_odds
std::vector< unsigned int > face_partition_data