29#ifdef DEAL_II_WITH_TBB
30# include <tbb/blocked_range.h>
31# include <tbb/parallel_for.h>
33# ifndef DEAL_II_TBB_WITH_ONEAPI
34# include <tbb/task_scheduler_init.h>
64 namespace MatrixFreeFunctions
66#if defined(DEAL_II_WITH_TBB) && !defined(DEAL_II_TBB_WITH_ONEAPI)
76 ActualCellWork(MFWorkerInterface **worker_pointer,
77 const unsigned int partition,
78 const TaskInfo & task_info)
80 , worker_pointer(worker_pointer)
82 , task_info(task_info)
85 ActualCellWork(MFWorkerInterface &worker,
86 const unsigned int partition,
87 const TaskInfo & task_info)
89 , worker_pointer(nullptr)
91 , task_info(task_info)
97 MFWorkerInterface *used_worker =
98 worker !=
nullptr ? worker : *worker_pointer;
100 used_worker->cell(partition);
102 if (task_info.face_partition_data.empty() ==
false)
104 used_worker->face(partition);
105 used_worker->boundary(partition);
110 MFWorkerInterface * worker;
111 MFWorkerInterface **worker_pointer;
113 const TaskInfo & task_info;
116 class CellWork :
public tbb::task
119 CellWork(MFWorkerInterface &worker,
120 const unsigned int partition,
121 const TaskInfo & task_info,
122 const bool is_blocked)
125 , is_blocked(is_blocked)
133 if (is_blocked ==
true)
134 tbb::empty_task::spawn(*dummy);
138 tbb::empty_task *dummy;
142 const bool is_blocked;
147 class PartitionWork :
public tbb::task
150 PartitionWork(MFWorkerInterface &function_in,
151 const unsigned int partition_in,
152 const TaskInfo & task_info_in,
153 const bool is_blocked_in =
false)
155 , function(function_in)
157 , task_info(task_info_in)
158 , is_blocked(is_blocked_in)
164 tbb::empty_task *root =
165 new (tbb::task::allocate_root()) tbb::empty_task;
166 const unsigned int evens = task_info.partition_evens[
partition];
167 const unsigned int odds = task_info.partition_odds[
partition];
168 const unsigned int n_blocked_workers =
169 task_info.partition_n_blocked_workers[
partition];
170 const unsigned int n_workers =
171 task_info.partition_n_workers[
partition];
172 std::vector<CellWork *> worker(n_workers);
173 std::vector<CellWork *> blocked_worker(n_blocked_workers);
175 root->set_ref_count(evens + 1);
176 for (
unsigned int j = 0; j < evens; ++j)
178 worker[j] =
new (root->allocate_child())
180 task_info.partition_row_index[partition] + 2 * j,
185 worker[j]->set_ref_count(2);
186 blocked_worker[j - 1]->dummy =
187 new (worker[j]->allocate_child()) tbb::empty_task;
188 tbb::task::spawn(*blocked_worker[j - 1]);
191 worker[j]->set_ref_count(1);
194 blocked_worker[j] =
new (worker[j]->allocate_child())
196 task_info.partition_row_index[partition] + 2 * j +
205 worker[evens] =
new (worker[j]->allocate_child())
207 task_info.partition_row_index[partition] +
211 tbb::task::spawn(*worker[evens]);
215 tbb::empty_task *child =
216 new (worker[j]->allocate_child()) tbb::empty_task();
217 tbb::task::spawn(*child);
222 root->wait_for_all();
223 root->destroy(*root);
224 if (is_blocked ==
true)
225 tbb::empty_task::spawn(*dummy);
229 tbb::empty_task *dummy;
232 MFWorkerInterface &function;
234 const TaskInfo & task_info;
235 const bool is_blocked;
247 CellWork(MFWorkerInterface &worker_in,
248 const TaskInfo & task_info_in,
249 const unsigned int partition_in)
251 , task_info(task_info_in)
256 operator()(
const tbb::blocked_range<unsigned int> &r)
const
258 const unsigned int start_index =
259 task_info.cell_partition_data[
partition] +
260 task_info.block_size * r.begin();
261 const unsigned int end_index =
262 std::min(start_index + task_info.block_size * (r.end() - r.begin()),
263 task_info.cell_partition_data[partition + 1]);
264 worker.cell(std::make_pair(start_index, end_index));
266 if (task_info.face_partition_data.empty() ==
false)
273 MFWorkerInterface &worker;
274 const TaskInfo & task_info;
280 class PartitionWork :
public tbb::task
283 PartitionWork(MFWorkerInterface &worker_in,
284 const unsigned int partition_in,
285 const TaskInfo & task_info_in,
286 const bool is_blocked_in)
290 , task_info(task_info_in)
291 , is_blocked(is_blocked_in)
297 const unsigned int n_chunks =
298 (task_info.cell_partition_data[
partition + 1] -
299 task_info.cell_partition_data[
partition] + task_info.block_size -
301 task_info.block_size;
302 parallel_for(tbb::blocked_range<unsigned int>(0, n_chunks, 1),
303 CellWork(worker, task_info, partition));
304 if (is_blocked ==
true)
305 tbb::empty_task::spawn(*dummy);
309 tbb::empty_task *dummy;
312 MFWorkerInterface &worker;
314 const TaskInfo & task_info;
315 const bool is_blocked;
322 class MPICommunication :
public tbb::task
325 MPICommunication(MFWorkerInterface &worker_in,
const bool do_compress)
327 , do_compress(do_compress)
333 if (do_compress ==
false)
334 worker.vector_update_ghosts_finish();
336 worker.vector_compress_start();
341 MFWorkerInterface &worker;
342 const bool do_compress;
364#if defined(DEAL_II_WITH_TBB) && !defined(DEAL_II_TBB_WITH_ONEAPI)
371 tbb::empty_task *root =
372 new (tbb::task::allocate_root()) tbb::empty_task;
373 root->set_ref_count(
evens + 1);
374 std::vector<partition::PartitionWork *> worker(
n_workers);
375 std::vector<partition::PartitionWork *> blocked_worker(
377 MPICommunication *worker_compr =
378 new (root->allocate_child()) MPICommunication(funct,
true);
379 worker_compr->set_ref_count(1);
380 for (
unsigned int j = 0; j <
evens; ++j)
384 worker[j] =
new (root->allocate_child())
385 partition::PartitionWork(funct, 2 * j, *
this,
false);
386 worker[j]->set_ref_count(2);
387 blocked_worker[j - 1]->dummy =
388 new (worker[j]->allocate_child()) tbb::empty_task;
389 tbb::task::spawn(*blocked_worker[j - 1]);
393 worker[j] =
new (worker_compr->allocate_child())
394 partition::PartitionWork(funct, 2 * j, *
this,
false);
395 worker[j]->set_ref_count(2);
396 MPICommunication *worker_dist =
397 new (worker[j]->allocate_child())
398 MPICommunication(funct,
false);
399 tbb::task::spawn(*worker_dist);
403 blocked_worker[j] =
new (worker[j]->allocate_child())
404 partition::PartitionWork(funct, 2 * j + 1, *
this,
true);
410 worker[
evens] =
new (worker[j]->allocate_child())
411 partition::PartitionWork(funct,
415 tbb::task::spawn(*worker[
evens]);
419 tbb::empty_task *child =
420 new (worker[j]->allocate_child()) tbb::empty_task();
421 tbb::task::spawn(*child);
426 root->wait_for_all();
427 root->destroy(*root);
443 tbb::empty_task *root =
444 new (tbb::task::allocate_root()) tbb::empty_task;
445 root->set_ref_count(
evens + 1);
450 std::vector<color::PartitionWork *> worker(
n_workers);
451 std::vector<color::PartitionWork *> blocked_worker(
453 unsigned int worker_index = 0, slice_index = 0;
454 int spawn_index_child = -2;
455 MPICommunication *worker_compr =
456 new (root->allocate_child()) MPICommunication(funct,
true);
457 worker_compr->set_ref_count(1);
458 for (
unsigned int part = 0;
463 worker[worker_index] =
464 new (worker_compr->allocate_child())
465 color::PartitionWork(funct,
470 worker[worker_index] =
new (root->allocate_child())
471 color::PartitionWork(funct,
479 worker[worker_index]->set_ref_count(1);
481 worker[worker_index] =
482 new (worker[worker_index - 1]->allocate_child())
483 color::PartitionWork(funct,
488 worker[worker_index]->set_ref_count(2);
491 blocked_worker[(part - 1) / 2]->dummy =
492 new (worker[worker_index]->allocate_child())
495 if (spawn_index_child == -1)
496 tbb::task::spawn(*blocked_worker[(part - 1) / 2]);
499 Assert(spawn_index_child >= 0,
501 tbb::task::spawn(*worker[spawn_index_child]);
503 spawn_index_child = -2;
507 MPICommunication *worker_dist =
508 new (worker[worker_index]->allocate_child())
509 MPICommunication(funct,
false);
510 tbb::task::spawn(*worker_dist);
518 blocked_worker[part / 2] =
519 new (worker[worker_index - 1]->allocate_child())
520 color::PartitionWork(funct,
527 blocked_worker[part / 2]->set_ref_count(1);
528 worker[worker_index] =
new (
529 blocked_worker[part / 2]->allocate_child())
530 color::PartitionWork(funct,
538 spawn_index_child = -1;
547 worker[worker_index]->set_ref_count(1);
550 worker[worker_index] =
551 new (worker[worker_index - 1]->allocate_child())
552 color::PartitionWork(funct,
557 spawn_index_child = worker_index;
562 tbb::empty_task *
final =
563 new (worker[worker_index - 1]->allocate_child())
565 tbb::task::spawn(*
final);
566 spawn_index_child = worker_index - 1;
572 tbb::task::spawn(*worker[spawn_index_child]);
574 root->wait_for_all();
575 root->destroy(*root);
588 tbb::empty_task *root =
589 new (tbb::task::allocate_root()) tbb::empty_task;
590 root->set_ref_count(2);
591 color::PartitionWork *worker =
592 new (root->allocate_child())
593 color::PartitionWork(funct,
color, *
this,
false);
594 tbb::empty_task::spawn(*worker);
595 root->wait_for_all();
596 root->destroy(*root);
688 template <
typename StreamType>
691 const std::size_t data_length)
const
698 out << memory_c.
min <<
"/" << memory_c.
avg <<
"/" << memory_c.
max;
699 out <<
" MB" << std::endl;
723 std::vector<unsigned int> &boundary_cells)
727 unsigned int fillup_needed =
735 std::vector<unsigned int> new_boundary_cells;
736 new_boundary_cells.reserve(boundary_cells.size());
738 unsigned int next_free_slot = 0, bound_index = 0;
739 while (fillup_needed > 0 && bound_index < boundary_cells.size())
741 if (next_free_slot < boundary_cells[bound_index])
745 if (next_free_slot + fillup_needed <=
746 boundary_cells[bound_index])
748 for (
unsigned int j =
749 boundary_cells[bound_index] - fillup_needed;
750 j < boundary_cells[bound_index];
752 new_boundary_cells.push_back(j);
759 for (
unsigned int j = next_free_slot;
760 j < boundary_cells[bound_index];
762 new_boundary_cells.push_back(j);
764 boundary_cells[bound_index] - next_free_slot;
767 new_boundary_cells.push_back(boundary_cells[bound_index]);
768 next_free_slot = boundary_cells[bound_index] + 1;
771 while (fillup_needed > 0 &&
772 (new_boundary_cells.size() == 0 ||
774 new_boundary_cells.push_back(new_boundary_cells.back() + 1);
775 while (bound_index < boundary_cells.size())
776 new_boundary_cells.push_back(boundary_cells[bound_index++]);
778 boundary_cells.swap(new_boundary_cells);
782 std::sort(boundary_cells.begin(), boundary_cells.end());
795 const std::vector<unsigned int> &cells_with_comm,
796 const unsigned int dofs_per_cell,
797 const bool categories_are_hp,
798 const std::vector<unsigned int> &cell_vectorization_categories,
799 const bool cell_vectorization_categories_strict,
800 const std::vector<unsigned int> &parent_relation,
801 std::vector<unsigned int> & renumbering,
802 std::vector<unsigned char> & incompletely_filled_vectorization)
828 unsigned int vectorization_length_bits = 0;
830 while ((my_length >>= 1) != 0u)
831 ++vectorization_length_bits;
832 const unsigned int n_lanes = 1 << vectorization_length_bits;
837 unsigned int n_categories = 1;
839 if (cell_vectorization_categories.empty() ==
false)
844 std::set<unsigned int> used_categories;
846 used_categories.insert(cell_vectorization_categories[i]);
847 std::vector<unsigned int> used_categories_vector(
848 used_categories.size());
850 for (
const auto &it : used_categories)
851 used_categories_vector[n_categories++] = it;
854 const unsigned int index =
855 std::lower_bound(used_categories_vector.begin(),
856 used_categories_vector.end(),
857 cell_vectorization_categories[i]) -
858 used_categories_vector.begin();
860 tight_category_map[i] =
index;
867 std::vector<std::vector<unsigned int>> renumbering_category(n_categories);
869 renumbering_category[tight_category_map[i]].push_back(i);
871 if (cell_vectorization_categories_strict ==
false && n_categories > 1)
872 for (
unsigned int j = n_categories - 1; j > 0; --j)
874 unsigned int lower_index = j - 1;
875 while ((renumbering_category[j].size() % n_lanes) != 0u)
877 while (((renumbering_category[j].size() % n_lanes) != 0u) &&
878 !renumbering_category[lower_index].empty())
880 renumbering_category[j].push_back(
881 renumbering_category[lower_index].back());
882 renumbering_category[lower_index].pop_back();
884 if (lower_index == 0)
895 std::vector<unsigned int> temporary_numbering;
897 (n_lanes - 1) * n_categories);
898 const unsigned int n_cells_per_parent =
899 std::count(parent_relation.begin(), parent_relation.end(), 0);
900 std::vector<unsigned int> category_size;
901 for (
unsigned int j = 0; j < n_categories; ++j)
903 std::vector<std::pair<unsigned int, unsigned int>> grouped_cells;
904 std::vector<unsigned int> other_cells;
905 for (
const unsigned int cell : renumbering_category[j])
906 if (parent_relation.empty() ||
908 other_cells.push_back(cell);
910 grouped_cells.emplace_back(parent_relation[cell], cell);
913 std::sort(grouped_cells.begin(), grouped_cells.end());
914 std::vector<unsigned int> n_cells_per_group;
915 unsigned int length = 0;
916 for (
unsigned int i = 0; i < grouped_cells.size(); ++i, ++length)
917 if (i > 0 && grouped_cells[i].
first != grouped_cells[i - 1].
first)
919 n_cells_per_group.push_back(length);
923 n_cells_per_group.push_back(length);
928 auto group_it = grouped_cells.begin();
929 for (
unsigned int length : n_cells_per_group)
930 if (length < n_cells_per_parent)
931 for (
unsigned int j = 0; j < length; ++j)
932 other_cells.push_back((group_it++)->second);
938 for (
unsigned int j = 0; j < length; ++j)
939 temporary_numbering.push_back((group_it++)->second);
943 std::sort(other_cells.begin(), other_cells.end());
944 temporary_numbering.insert(temporary_numbering.end(),
948 while (temporary_numbering.size() % n_lanes != 0)
951 category_size.push_back(temporary_numbering.size());
955 std::vector<bool> batch_with_comm(temporary_numbering.size() / n_lanes,
957 std::vector<unsigned int> temporary_numbering_inverse(
n_active_cells);
958 for (
unsigned int i = 0; i < temporary_numbering.size(); ++i)
960 temporary_numbering_inverse[temporary_numbering[i]] = i;
961 for (
const unsigned int cell : cells_with_comm)
962 batch_with_comm[temporary_numbering_inverse[cell] / n_lanes] =
true;
968 std::vector<std::array<unsigned int, 3>> batch_order;
969 std::vector<std::array<unsigned int, 3>> batch_order_comm;
970 for (
unsigned int i = 0; i < temporary_numbering.size(); i += n_lanes)
972 unsigned int max_index = 0;
973 for (
unsigned int j = 0; j < n_lanes; ++j)
975 max_index =
std::max(temporary_numbering[i + j], max_index);
976 const unsigned int category_hp =
978 std::upper_bound(category_size.begin(), category_size.end(), i) -
979 category_size.begin() :
981 const std::array<unsigned int, 3> next{{category_hp, max_index, i}};
982 if (batch_with_comm[i / n_lanes])
983 batch_order_comm.emplace_back(next);
985 batch_order.emplace_back(next);
988 std::sort(batch_order.begin(), batch_order.end());
989 std::sort(batch_order_comm.begin(), batch_order_comm.end());
996 std::vector<unsigned int> blocks;
999 if (batch_order.empty())
1000 std::swap(batch_order_comm, batch_order);
1003 blocks = {0,
static_cast<unsigned int>(batch_order.size())};
1008 const unsigned int comm_begin = batch_order.size() / 2;
1009 batch_order.insert(batch_order.begin() + comm_begin,
1010 batch_order_comm.begin(),
1011 batch_order_comm.end());
1012 const unsigned int comm_end = comm_begin + batch_order_comm.size();
1013 const unsigned int end = batch_order.size();
1014 blocks = {0, comm_begin, comm_end, end};
1018 const unsigned int n_cell_batches = batch_order.size();
1019 const unsigned int n_ghost_batches =
1021 incompletely_filled_vectorization.resize(n_cell_batches +
1027 renumbering.clear();
1031 unsigned int counter = 0;
1032 for (
unsigned int block = 0; block < blocks.size() - 1; ++block)
1034 const unsigned int grain_size =
1035 std::max((2048U / dofs_per_cell) / 8 * 4, 2U);
1036 for (
unsigned int k = blocks[block]; k < blocks[block + 1];
1039 std::min(k + grain_size, blocks[block + 1]));
1043 for (
unsigned int k = blocks[block]; k < blocks[block + 1]; ++k)
1045 const unsigned int pos = batch_order[k][2];
1047 for (; j < n_lanes && temporary_numbering[pos + j] !=
1050 renumbering[counter++] = temporary_numbering[pos + j];
1052 incompletely_filled_vectorization[k] = j;
1062 if (!cell_vectorization_categories.empty())
1065 renumbering[cell] = cell;
1068 incompletely_filled_vectorization.back() =
n_ghost_cells % n_lanes;
1073 std::vector<unsigned int> renumber_cpy(renumbering);
1074 std::sort(renumber_cpy.begin(), renumber_cpy.end());
1075 for (
unsigned int i = 0; i < renumber_cpy.size(); ++i)
1084 const std::vector<unsigned int> &boundary_cells,
1085 std::vector<unsigned int> & renumbering,
1086 std::vector<unsigned char> & incompletely_filled_vectorization)
1088 const unsigned int n_cell_batches =
1090 const unsigned int n_ghost_slots =
1092 incompletely_filled_vectorization.resize(n_cell_batches + n_ghost_slots);
1094 incompletely_filled_vectorization[n_cell_batches - 1] =
1098 incompletely_filled_vectorization[n_cell_batches + n_ghost_slots - 1] =
1102 std::vector<unsigned int> reverse_numbering(
1104 for (
unsigned int j = 0; j < boundary_cells.size(); ++j)
1105 reverse_numbering[boundary_cells[j]] = j;
1106 unsigned int counter = boundary_cells.size();
1109 reverse_numbering[j] = counter++;
1116 renumbering.push_back(j);
1124 const unsigned int n_macro_boundary_cells =
1128 (n_cell_batches - n_macro_boundary_cells) / 2);
1130 n_macro_boundary_cells);
1161 const unsigned int minimum_parallel_grain_size = 200;
1162 if (dofs_per_cell *
block_size < minimum_parallel_grain_size)
1163 block_size = (minimum_parallel_grain_size / dofs_per_cell + 1);
1168 1 <<
static_cast<unsigned int>(std::log2(
block_size + 1));
1179 std::vector<unsigned int> & renumbering,
1180 std::vector<unsigned char> &irregular_cells,
1184 if (n_cell_batches == 0)
1189 unsigned int partition = 0, counter = 0;
1203 std::vector<unsigned int> cell_partition(
n_blocks,
1208 std::vector<unsigned int> partition_list(
n_blocks, 0);
1209 std::vector<unsigned int> partition_color_list(
n_blocks, 0);
1212 std::vector<unsigned int> partition_size(2, 0);
1218 unsigned int cluster_size = 1;
1234 partition_color_list);
1236 partition_list = renumbering;
1241 std::vector<unsigned int> sorted_pc_list(partition_color_list);
1242 std::sort(sorted_pc_list.begin(), sorted_pc_list.end());
1243 for (
unsigned int i = 0; i < sorted_pc_list.size(); ++i)
1250 std::vector<unsigned int> block_start(n_cell_batches + 1);
1251 std::vector<unsigned char> irregular(n_cell_batches);
1253 unsigned int mcell_start = 0;
1255 for (
unsigned int block = 0; block <
n_blocks; ++block)
1257 block_start[block + 1] = block_start[block];
1258 for (
unsigned int mcell = mcell_start;
1262 unsigned int n_comp = (irregular_cells[mcell] > 0) ?
1263 irregular_cells[mcell] :
1265 block_start[block + 1] += n_comp;
1271 unsigned int counter_macro = 0;
1272 unsigned int block_size_last =
1274 if (block_size_last == 0)
1277 unsigned int tick = 0;
1278 for (
unsigned int block = 0; block <
n_blocks; ++block)
1280 unsigned int present_block = partition_color_list[block];
1281 for (
unsigned int cell = block_start[present_block];
1282 cell < block_start[present_block + 1];
1284 renumbering[counter++] = partition_list[cell];
1285 unsigned int this_block_size =
1293 for (
unsigned int j = 0; j < this_block_size; ++j)
1294 irregular[counter_macro++] =
1295 irregular_cells[present_block *
block_size + j];
1300 irregular_cells.swap(irregular);
1307 std::vector<unsigned int> sorted_renumbering(renumbering);
1308 std::sort(sorted_renumbering.begin(), sorted_renumbering.end());
1309 for (
unsigned int i = 0; i < sorted_renumbering.size(); ++i)
1325 const std::vector<unsigned int> &cell_active_fe_index,
1327 std::vector<unsigned int> & renumbering,
1328 std::vector<unsigned char> & irregular_cells,
1332 if (n_cell_batches == 0)
1342 connectivity_blocks);
1354 std::vector<unsigned int> cell_partition(
n_blocks,
1360 std::vector<unsigned int> partition_list(
n_blocks, 0);
1361 std::vector<unsigned int> partition_2layers_list(
n_blocks, 0);
1364 std::vector<unsigned int> partition_size(2, 0);
1366 unsigned int partition = 0;
1372 unsigned int cluster_size = 1;
1399 cell_active_fe_index,
1406 partition_2layers_list,
1416 partition_2layers_list);
1422 std::vector<unsigned int> sorted_pc_list(partition_2layers_list);
1423 std::sort(sorted_pc_list.begin(), sorted_pc_list.end());
1424 for (
unsigned int i = 0; i < sorted_pc_list.size(); ++i)
1431 renumbering_in.swap(renumbering);
1437 for (
unsigned int j = 0; j < renumbering.size(); ++j)
1438 renumbering[j] = renumbering_in[partition_2layers_list[j]];
1447 std::vector<unsigned int> block_start(n_cell_batches + 1);
1448 std::vector<unsigned char> irregular(n_cell_batches);
1450 unsigned int counter = 0;
1451 unsigned int mcell_start = 0;
1453 for (
unsigned int block = 0; block <
n_blocks; ++block)
1455 block_start[block + 1] = block_start[block];
1456 for (
unsigned int mcell = mcell_start;
1460 unsigned int n_comp = (irregular_cells[mcell] > 0) ?
1461 irregular_cells[mcell] :
1463 block_start[block + 1] += n_comp;
1469 unsigned int counter_macro = 0;
1470 unsigned int block_size_last =
1472 if (block_size_last == 0)
1475 unsigned int tick = 0;
1476 for (
unsigned int block = 0; block <
n_blocks; ++block)
1478 unsigned int present_block = partition_2layers_list[block];
1479 for (
unsigned int cell = block_start[present_block];
1480 cell < block_start[present_block + 1];
1482 renumbering[counter++] = renumbering_in[cell];
1483 unsigned int this_block_size =
1491 for (
unsigned int j = 0; j < this_block_size; ++j)
1492 irregular[counter_macro++] =
1493 irregular_cells[present_block *
block_size + j];
1498 irregular_cells.swap(irregular);
1504 std::vector<unsigned int> sorted_renumbering(renumbering);
1505 std::sort(sorted_renumbering.begin(), sorted_renumbering.end());
1506 for (
unsigned int i = 0; i < sorted_renumbering.size(); ++i)
1520 const std::vector<unsigned int> &cell_active_fe_index,
1522 std::vector<unsigned int> & renumbering,
1523 std::vector<unsigned char> & irregular_cells,
1527 if (n_cell_batches == 0)
1544 std::vector<unsigned int> partition_partition_list(
n_active_cells, 0);
1547 std::vector<unsigned int> partition_size(2, 0);
1549 unsigned int partition = 0;
1563 cell_active_fe_index,
1570 partition_partition_list,
1573 partition_list.swap(renumbering);
1575 for (
unsigned int j = 0; j < renumbering.size(); ++j)
1576 renumbering[j] = partition_list[partition_partition_list[j]];
1588 const std::vector<unsigned char> &irregular_cells,
1592 std::vector<std::vector<unsigned int>> cell_blocks(
n_blocks);
1594 unsigned int cell = 0;
1595 for (
unsigned int i = 0, mcell = 0; i <
n_blocks; ++i)
1597 for (
unsigned int c = 0;
1601 unsigned int ncomp = (irregular_cells[mcell] > 0) ?
1602 irregular_cells[mcell] :
1604 for (
unsigned int c = 0; c < ncomp; ++c, ++cell)
1606 cell_blocks[i].push_back(cell);
1607 touched_cells[cell] = i;
1612 for (
unsigned int i = 0; i < cell_blocks.size(); ++i)
1613 for (
unsigned int col = 0; col < cell_blocks[i].size(); ++col)
1616 connectivity_cells.
begin(cell_blocks[i][col]);
1617 it != connectivity_cells.
end(cell_blocks[i][col]);
1620 if (touched_cells[it->column()] != i)
1621 connectivity_blocks.
add(i, touched_cells[it->column()]);
1633 const std::vector<unsigned int> &cell_active_fe_index,
1634 const unsigned int partition,
1635 const unsigned int cluster_size,
1637 const std::vector<unsigned int> &cell_partition,
1638 const std::vector<unsigned int> &partition_list,
1639 const std::vector<unsigned int> &partition_size,
1640 std::vector<unsigned int> & partition_partition_list,
1641 std::vector<unsigned char> & irregular_cells)
1644 const unsigned int n_ghost_slots =
1648 std::vector<unsigned int> neighbor_list;
1651 std::vector<unsigned int> neighbor_neighbor_list;
1655 irregular_cells.back() = 0;
1658 unsigned int max_fe_index = 0;
1659 for (
const unsigned int fe_index : cell_active_fe_index)
1660 max_fe_index =
std::max(fe_index, max_fe_index);
1666 unsigned int n_cell_batches_before = 0;
1672 std::vector<unsigned int> cell_partition_l2(
1678 unsigned int counter = 0;
1679 unsigned int missing_macros;
1680 for (
unsigned int part = 0; part < partition; ++part)
1682 neighbor_neighbor_list.resize(0);
1683 neighbor_list.resize(0);
1685 unsigned int partition_l2 = 0;
1686 unsigned int start_up = partition_size[part];
1687 unsigned int partition_counter = 0;
1690 if (neighbor_list.size() == 0)
1693 partition_counter = 0;
1694 for (
unsigned int j = start_up;
1695 j < partition_size[part + 1];
1697 if (cell_partition[partition_list[j]] == part &&
1698 cell_partition_l2[partition_list[j]] ==
1703 partition_counter = 1;
1707 cell_partition_l2[partition_list[start_up]] =
1709 neighbor_neighbor_list.push_back(
1710 partition_list[start_up]);
1711 partition_partition_list[counter++] =
1712 partition_list[start_up];
1719 partition_counter = 0;
1720 for (
const unsigned int neighbor : neighbor_list)
1722 Assert(cell_partition[neighbor] == part,
1724 Assert(cell_partition_l2[neighbor] == partition_l2 - 1,
1726 auto neighbor_it = connectivity.
begin(neighbor);
1727 const auto end_it = connectivity.
end(neighbor);
1728 for (; neighbor_it != end_it; ++neighbor_it)
1730 if (cell_partition[neighbor_it->column()] == part &&
1731 cell_partition_l2[neighbor_it->column()] ==
1734 cell_partition_l2[neighbor_it->column()] =
1736 neighbor_neighbor_list.push_back(
1737 neighbor_it->column());
1738 partition_partition_list[counter++] =
1739 neighbor_it->column();
1740 partition_counter++;
1745 if (partition_counter > 0)
1747 int index_before = neighbor_neighbor_list.size(),
1748 index = index_before;
1753 std::vector<unsigned int> remaining_per_cell_batch(
1755 std::vector<std::vector<unsigned int>>
1756 renumbering_fe_index;
1759 if (hp_bool ==
true)
1761 renumbering_fe_index.resize(max_fe_index + 1);
1762 for (cell = counter - partition_counter;
1766 renumbering_fe_index
1767 [cell_active_fe_index.empty() ?
1769 cell_active_fe_index
1770 [partition_partition_list[cell]]]
1771 .push_back(partition_partition_list[cell]);
1774 for (
unsigned int j = 0; j < max_fe_index + 1; ++j)
1776 remaining_per_cell_batch[j] =
1777 renumbering_fe_index[j].size() %
1779 if (remaining_per_cell_batch[j] != 0)
1782 ((renumbering_fe_index[j].size() +
1789 remaining_per_cell_batch.resize(1);
1790 remaining_per_cell_batch[0] =
1794 if (remaining_per_cell_batch[0] != 0)
1801 cluster_size - (missing_macros % cluster_size);
1804 while (missing_macros > 0 || filled ==
false)
1808 index = neighbor_neighbor_list.size();
1809 if (
index == index_before)
1811 if (missing_macros != 0)
1813 neighbor_neighbor_list.resize(0);
1818 index_before =
index;
1821 unsigned int additional =
1822 neighbor_neighbor_list[
index];
1833 for (; neighbor != end; ++neighbor)
1835 if (cell_partition[neighbor->
column()] == part &&
1836 cell_partition_l2[neighbor->
column()] ==
1839 unsigned int this_index = 0;
1840 if (hp_bool ==
true)
1842 cell_active_fe_index.empty() ?
1844 cell_active_fe_index[neighbor
1851 if (missing_macros > 0 ||
1852 remaining_per_cell_batch[this_index] > 0)
1854 cell_partition_l2[neighbor->
column()] =
1856 neighbor_neighbor_list.push_back(
1858 if (hp_bool ==
true)
1859 renumbering_fe_index[this_index]
1860 .push_back(neighbor->
column());
1861 partition_partition_list[counter] =
1864 partition_counter++;
1865 if (remaining_per_cell_batch
1866 [this_index] == 0 &&
1869 remaining_per_cell_batch[this_index]++;
1870 if (remaining_per_cell_batch
1874 remaining_per_cell_batch[this_index] =
1877 if (missing_macros == 0)
1880 for (
unsigned int fe_ind = 0;
1881 fe_ind < max_fe_index + 1;
1883 if (remaining_per_cell_batch
1893 if (hp_bool ==
true)
1898 cell = counter - partition_counter;
1899 for (
unsigned int j = 0; j < max_fe_index + 1; ++j)
1901 for (
const unsigned int jj :
1902 renumbering_fe_index[j])
1903 renumbering[cell++] = jj;
1904 if (renumbering_fe_index[j].size() %
1907 irregular_cells[renumbering_fe_index[j].size() /
1909 n_cell_batches_before] =
1910 renumbering_fe_index[j].size() %
1912 n_cell_batches_before +=
1913 (renumbering_fe_index[j].size() +
1916 renumbering_fe_index[j].resize(0);
1921 n_cell_batches_before +=
1925 irregular_cells[n_cell_batches_before] =
1927 n_cell_batches_before++;
1934 neighbor_list = neighbor_neighbor_list;
1935 neighbor_neighbor_list.resize(0);
1941 if (hp_bool ==
true)
1943 partition_partition_list.swap(renumbering);
1954 const unsigned int partition,
1955 const std::vector<unsigned int> &cell_partition,
1956 const std::vector<unsigned int> &partition_list,
1957 const std::vector<unsigned int> &partition_size,
1958 std::vector<unsigned int> & partition_color_list)
1961 std::vector<unsigned int> cell_color(
n_blocks, n_cell_batches);
1962 std::vector<bool> color_finder;
1966 unsigned int color_counter = 0, index_counter = 0;
1967 for (
unsigned int part = 0; part < partition; ++part)
1970 unsigned int max_color = 0;
1971 for (
unsigned int k = partition_size[part];
1972 k < partition_size[part + 1];
1975 unsigned int cell = partition_list[k];
1976 unsigned int n_neighbors = connectivity.
row_length(cell);
1980 color_finder.resize(n_neighbors + 1);
1981 for (
unsigned int j = 0; j <= n_neighbors; ++j)
1982 color_finder[j] =
true;
1984 connectivity.
begin(cell),
1985 end = connectivity.
end(cell);
1986 for (; neighbor != end; ++neighbor)
1990 if (cell_partition[neighbor->
column()] == part &&
1991 cell_color[neighbor->
column()] <= n_neighbors)
1992 color_finder[cell_color[neighbor->
column()]] =
false;
1995 cell_color[cell] = 0;
1996 while (color_finder[cell_color[cell]] ==
false)
1998 if (cell_color[cell] > max_color)
1999 max_color = cell_color[cell];
2008 for (
unsigned int k = partition_size[part];
2009 k < partition_size[part + 1];
2012 unsigned int cell = partition_list[k];
2013 if (cell_color[cell] ==
color)
2015 partition_color_list[color_counter++] = cell;
2029 const unsigned int cluster_size,
2030 std::vector<unsigned int> & cell_partition,
2031 std::vector<unsigned int> & partition_list,
2032 std::vector<unsigned int> & partition_size,
2033 unsigned int & partition)
const
2042 std::vector<unsigned int> neighbor_list;
2045 std::vector<unsigned int> neighbor_neighbor_list;
2055 unsigned int counter = 0;
2056 unsigned int start_nonboundary =
2063 if (n_cell_batches == 0)
2066 start_nonboundary = n_cell_batches;
2081 unsigned int start_up = 0;
2083 unsigned int remainder = cluster_size;
2091 if (start_nonboundary > 0)
2093 for (
unsigned int cell = 0; cell < start_nonboundary; ++cell)
2095 const unsigned int cell_nn = cell;
2096 cell_partition[cell_nn] = partition;
2097 neighbor_list.push_back(cell_nn);
2098 partition_list[counter++] = cell_nn;
2099 partition_size.back()++;
2101 start_nonboundary = 0;
2102 remainder -= (start_nonboundary % cluster_size);
2103 if (remainder == cluster_size)
2110 cell_partition[start_up] = partition;
2111 neighbor_list.push_back(start_up);
2112 partition_list[counter++] = start_up;
2113 partition_size.back()++;
2116 if (remainder == cluster_size)
2119 int index_before = neighbor_list.size(),
index = index_before,
2121 while (remainder > 0)
2123 if (
index == index_stop)
2125 index = neighbor_list.size();
2126 if (
index == index_before)
2128 neighbor_list.resize(0);
2131 index_stop = index_before;
2132 index_before =
index;
2135 unsigned int additional = neighbor_list[
index];
2137 connectivity.
begin(additional),
2139 connectivity.
end(additional);
2140 for (; neighbor != end; ++neighbor)
2142 if (cell_partition[neighbor->
column()] ==
2145 partition_size.back()++;
2146 cell_partition[neighbor->
column()] = partition;
2147 neighbor_list.push_back(neighbor->
column());
2148 partition_list[counter++] = neighbor->
column();
2156 while (neighbor_list.size() > 0)
2161 unsigned int partition_counter = 0;
2164 partition_size.push_back(partition_size.back());
2168 for (
const unsigned int cell : neighbor_list)
2170 Assert(cell_partition[cell] == partition - 1,
2172 auto neighbor = connectivity.
begin(cell);
2173 const auto end = connectivity.
end(cell);
2174 for (; neighbor != end; ++neighbor)
2176 if (cell_partition[neighbor->column()] ==
2179 partition_size.back()++;
2180 cell_partition[neighbor->column()] = partition;
2184 neighbor_neighbor_list.push_back(neighbor->column());
2185 partition_list[counter++] = neighbor->column();
2186 partition_counter++;
2190 remainder = cluster_size - (partition_counter % cluster_size);
2191 if (remainder == cluster_size)
2194 int index_before = neighbor_neighbor_list.size(),
2195 index = index_before;
2196 while (remainder > 0)
2198 if (
index == index_stop)
2200 index = neighbor_neighbor_list.size();
2201 if (
index == index_before)
2203 neighbor_neighbor_list.resize(0);
2206 index_stop = index_before;
2207 index_before =
index;
2210 unsigned int additional = neighbor_neighbor_list[
index];
2214 end = connectivity.
end(
2216 for (; neighbor != end; ++neighbor)
2218 if (cell_partition[neighbor->
column()] ==
2221 partition_size.back()++;
2222 cell_partition[neighbor->
column()] = partition;
2223 neighbor_neighbor_list.push_back(neighbor->
column());
2224 partition_list[counter++] = neighbor->
column();
2232 neighbor_list = neighbor_neighbor_list;
2233 neighbor_neighbor_list.resize(0);
2239 for (
unsigned int j = start_up; j <
n_blocks; ++j)
2245 remainder = cluster_size;
2259 evens = (partition + 1) / 2;
2260 odds = partition / 2;
2268 for (
unsigned int part = 0; part < partition; ++part)
2289internal::MatrixFreeFunctions::TaskInfo::print_memory_statistics<std::ostream>(
2291 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)
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
std::vector< Integer > invert_permutation(const std::vector< Integer > &permutation)
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