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> 24 #include <deal.II/matrix_free/task_info.h> 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> 37 DEAL_II_NAMESPACE_OPEN
44 namespace MatrixFreeFunctions
46 #ifdef DEAL_II_WITH_THREADS 56 ActualCellWork(MFWorkerInterface **worker_pointer,
57 const unsigned int partition,
58 const TaskInfo & task_info)
60 , worker_pointer(worker_pointer)
62 , task_info(task_info)
65 ActualCellWork(MFWorkerInterface &worker,
66 const unsigned int partition,
67 const TaskInfo & task_info)
69 , worker_pointer(nullptr)
71 , task_info(task_info)
77 MFWorkerInterface *used_worker =
78 worker !=
nullptr ? worker : *worker_pointer;
81 std::make_pair(task_info.cell_partition_data[partition],
82 task_info.cell_partition_data[partition + 1]));
84 if (task_info.face_partition_data.empty() ==
false)
87 std::make_pair(task_info.face_partition_data[partition],
88 task_info.face_partition_data[partition + 1]));
90 used_worker->boundary(std::make_pair(
91 task_info.boundary_partition_data[partition],
92 task_info.boundary_partition_data[partition + 1]));
97 MFWorkerInterface * worker;
98 MFWorkerInterface **worker_pointer;
100 const TaskInfo & task_info;
103 class CellWork :
public tbb::task
106 CellWork(MFWorkerInterface &worker,
107 const unsigned int partition,
108 const TaskInfo & task_info,
109 const bool is_blocked)
112 , is_blocked(is_blocked)
120 if (is_blocked ==
true)
121 tbb::empty_task::spawn(*dummy);
125 tbb::empty_task *dummy;
129 const bool is_blocked;
134 class PartitionWork :
public tbb::task
137 PartitionWork(MFWorkerInterface &function_in,
138 const unsigned int partition_in,
139 const TaskInfo & task_info_in,
140 const bool is_blocked_in =
false)
142 , function(function_in)
144 , task_info(task_info_in)
145 , is_blocked(is_blocked_in)
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);
162 root->set_ref_count(evens + 1);
163 for (
unsigned int j = 0; j < evens; j++)
165 worker[j] =
new (root->allocate_child())
167 task_info.partition_row_index[partition] + 2 * j,
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]);
178 worker[j]->set_ref_count(1);
181 blocked_worker[j] =
new (worker[j]->allocate_child())
183 task_info.partition_row_index[partition] + 2 * j +
192 worker[evens] =
new (worker[j]->allocate_child())
194 task_info.partition_row_index[partition] +
198 tbb::task::spawn(*worker[evens]);
202 tbb::empty_task *child =
203 new (worker[j]->allocate_child()) tbb::empty_task();
204 tbb::task::spawn(*child);
209 root->wait_for_all();
210 root->destroy(*root);
211 if (is_blocked ==
true)
212 tbb::empty_task::spawn(*dummy);
216 tbb::empty_task *dummy;
219 MFWorkerInterface &
function;
221 const TaskInfo & task_info;
222 const bool is_blocked;
234 CellWork(MFWorkerInterface &worker_in,
235 const TaskInfo & task_info_in,
236 const unsigned int partition_in)
238 , task_info(task_info_in)
243 operator()(
const tbb::blocked_range<unsigned int> &r)
const 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));
253 if (task_info.face_partition_data.empty() ==
false)
260 MFWorkerInterface &worker;
261 const TaskInfo & task_info;
267 class PartitionWork :
public tbb::task
270 PartitionWork(MFWorkerInterface &worker_in,
271 const unsigned int partition_in,
272 const TaskInfo & task_info_in,
273 const bool is_blocked_in)
277 , task_info(task_info_in)
278 , is_blocked(is_blocked_in)
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 -
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);
296 tbb::empty_task *dummy;
299 MFWorkerInterface &worker;
301 const TaskInfo & task_info;
302 const bool is_blocked;
309 class MPICommunication :
public tbb::task
312 MPICommunication(MFWorkerInterface &worker_in,
const bool do_compress)
314 , do_compress(do_compress)
320 if (do_compress ==
false)
321 worker.vector_update_ghosts_finish();
323 worker.vector_compress_start();
328 MFWorkerInterface &worker;
329 const bool do_compress;
332 #endif // DEAL_II_WITH_THREADS 341 #ifdef DEAL_II_WITH_THREADS 346 if (
scheme == partition_partition)
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++)
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]);
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);
380 blocked_worker[j] =
new (worker[j]->allocate_child())
381 partition::PartitionWork(funct, 2 * j + 1, *
this,
true);
387 worker[
evens] =
new (worker[j]->allocate_child())
388 partition::PartitionWork(funct,
392 tbb::task::spawn(*worker[
evens]);
396 tbb::empty_task *child =
397 new (worker[j]->allocate_child()) tbb::empty_task();
398 tbb::task::spawn(*child);
403 root->wait_for_all();
404 root->destroy(*root);
412 tbb::empty_task *root =
413 new (tbb::task::allocate_root()) tbb::empty_task;
414 root->set_ref_count(
evens + 1);
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;
432 worker[worker_index] =
433 new (worker_compr->allocate_child())
434 color::PartitionWork(funct,
439 worker[worker_index] =
new (root->allocate_child())
440 color::PartitionWork(funct,
448 worker[worker_index]->set_ref_count(1);
450 worker[worker_index] =
451 new (worker[worker_index - 1]->allocate_child())
452 color::PartitionWork(funct,
457 worker[worker_index]->set_ref_count(2);
460 blocked_worker[(part - 1) / 2]->dummy =
461 new (worker[worker_index]->allocate_child())
464 if (spawn_index_child == -1)
465 tbb::task::spawn(*blocked_worker[(part - 1) / 2]);
468 Assert(spawn_index_child >= 0,
470 tbb::task::spawn(*worker[spawn_index_child]);
472 spawn_index_child = -2;
476 MPICommunication *worker_dist =
477 new (worker[worker_index]->allocate_child())
478 MPICommunication(funct,
false);
479 tbb::task::spawn(*worker_dist);
487 blocked_worker[part / 2] =
488 new (worker[worker_index - 1]->allocate_child())
489 color::PartitionWork(funct,
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,
507 spawn_index_child = -1;
516 worker[worker_index]->set_ref_count(1);
519 worker[worker_index] =
520 new (worker[worker_index - 1]->allocate_child())
521 color::PartitionWork(funct,
526 spawn_index_child = worker_index;
531 tbb::empty_task *
final =
532 new (worker[worker_index - 1]->allocate_child())
534 tbb::task::spawn(*
final);
535 spawn_index_child = worker_index - 1;
541 tbb::task::spawn(*worker[spawn_index_child]);
543 root->wait_for_all();
544 root->destroy(*root);
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);
653 template <
typename StreamType>
656 const std::size_t data_length)
const 663 out << memory_c.
min <<
"/" << memory_c.
avg <<
"/" << memory_c.
max;
664 out <<
" MB" << std::endl;
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)
699 unsigned int fillup_needed =
707 std::vector<unsigned int> new_boundary_cells;
708 new_boundary_cells.reserve(boundary_cells.size());
710 unsigned int next_free_slot = 0, bound_index = 0;
711 while (fillup_needed > 0 && bound_index < boundary_cells.size())
713 if (next_free_slot < boundary_cells[bound_index])
717 if (next_free_slot + fillup_needed <=
718 boundary_cells[bound_index])
720 for (
unsigned int j =
721 boundary_cells[bound_index] - fillup_needed;
722 j < boundary_cells[bound_index];
724 new_boundary_cells.push_back(j);
731 for (
unsigned int j = next_free_slot;
732 j < boundary_cells[bound_index];
734 new_boundary_cells.push_back(j);
736 boundary_cells[bound_index] - next_free_slot;
739 new_boundary_cells.push_back(boundary_cells[bound_index]);
740 next_free_slot = boundary_cells[bound_index] + 1;
743 while (fillup_needed > 0 &&
744 (new_boundary_cells.size() == 0 ||
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++]);
750 boundary_cells.swap(new_boundary_cells);
754 std::sort(boundary_cells.begin(), boundary_cells.end());
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)
775 const unsigned int n_macro_cells =
777 const unsigned int n_ghost_slots =
779 const unsigned int n_boundary_cells = boundary_cells.size();
781 incompletely_filled_vectorization.resize(n_macro_cells + n_ghost_slots);
796 for (
unsigned int i = 0; i < n_boundary_cells; ++i)
797 cell_marked[boundary_cells[i]] = 2;
803 const unsigned int n_second_slot =
806 unsigned int count = 0;
807 for (
const unsigned int cell : cells_close_to_boundary)
808 if (cell_marked[cell] == 0)
810 cell_marked[cell] = count < n_second_slot ? 1 : 3;
816 if (cell_marked[c] == 0)
822 if (cell_marked[c] == 0)
825 if (cell_marked[c] == 0)
829 std::fill(cell_marked.begin(), cell_marked.end(), 1);
831 for (
const unsigned char marker : cell_marked)
837 unsigned int n_categories = 1;
838 std::vector<unsigned int> tight_category_map;
839 if (cell_vectorization_categories.empty() ==
false)
848 std::set<unsigned int> used_categories;
850 used_categories.insert(cell_vectorization_categories[i]);
851 std::vector<unsigned int> used_categories_vector(
852 used_categories.size());
854 for (
const auto &it : used_categories)
855 used_categories_vector[n_categories++] = it;
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();
864 tight_category_map[i] = index;
868 incompletely_filled_vectorization.resize(
869 incompletely_filled_vectorization.size() + 4 * n_categories);
871 else if (cells_close_to_boundary.empty())
877 for (
const unsigned int cell : cells_close_to_boundary)
878 tight_category_map[cell] = 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)
890 if (cell_marked[i] == block)
891 renumbering_category[tight_category_map[i]].push_back(i);
895 if (cell_vectorization_categories_strict ==
false && n_categories > 1)
896 for (
unsigned int j = n_categories - 1; j > 0; --j)
898 unsigned int lower_index = j - 1;
901 while (renumbering_category[j].size() %
903 !renumbering_category[lower_index].empty())
905 renumbering_category[j].push_back(
906 renumbering_category[lower_index].back());
907 renumbering_category[lower_index].pop_back();
909 if (lower_index == 0)
917 for (
unsigned int j = 0; j < n_categories; ++j)
919 for (
const unsigned int cell : renumbering_category[j])
920 renumbering[counter++] = cell;
921 unsigned int remainder =
924 incompletely_filled_vectorization
926 n_cells] = remainder;
927 const unsigned int n_my_macro_cells =
930 renumbering_category[j].clear();
935 std::max((2048U / dofs_per_cell) / 8 * 4, 2U);
937 for (
unsigned int k = 0; k < n_my_macro_cells; k +=
block_size)
939 n_cells + std::min(k +
block_size, n_my_macro_cells));
942 n_cells += n_my_macro_cells;
945 if (block == 3 || (block == 1 &&
n_procs == 1))
948 if (cell_vectorization_categories_strict ==
true)
966 const std::vector<unsigned int> &boundary_cells,
967 std::vector<unsigned int> & renumbering,
968 std::vector<unsigned char> & incompletely_filled_vectorization)
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);
976 incompletely_filled_vectorization[n_macro_cells - 1] =
980 incompletely_filled_vectorization[n_macro_cells + n_ghost_slots - 1] =
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();
991 reverse_numbering[j] = counter++;
998 renumbering.push_back(j);
1006 const unsigned int n_macro_boundary_cells =
1010 (n_macro_cells - n_macro_boundary_cells) / 2);
1012 n_macro_boundary_cells);
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);
1050 1 <<
static_cast<unsigned int>(std::log2(
block_size + 1));
1061 std::vector<unsigned int> & renumbering,
1062 std::vector<unsigned char> &irregular_cells,
1066 if (n_macro_cells == 0)
1071 unsigned int partition = 0, counter = 0;
1085 std::vector<unsigned int> cell_partition(
n_blocks,
1090 std::vector<unsigned int> partition_list(
n_blocks, 0);
1091 std::vector<unsigned int> partition_color_list(
n_blocks, 0);
1094 std::vector<unsigned int> partition_size(2, 0);
1100 unsigned int cluster_size = 1;
1116 partition_color_list);
1118 partition_list = renumbering;
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)
1132 std::vector<unsigned int> block_start(n_macro_cells + 1);
1133 std::vector<unsigned char> irregular(n_macro_cells);
1135 unsigned int mcell_start = 0;
1137 for (
unsigned int block = 0; block <
n_blocks; block++)
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);
1144 unsigned int n_comp = (irregular_cells[mcell] > 0) ?
1145 irregular_cells[mcell] :
1147 block_start[block + 1] += n_comp;
1153 unsigned int counter_macro = 0;
1154 unsigned int block_size_last =
1156 if (block_size_last == 0)
1159 unsigned int tick = 0;
1160 for (
unsigned int block = 0; block <
n_blocks; block++)
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];
1166 renumbering[counter++] = partition_list[cell];
1167 unsigned int this_block_size =
1175 for (
unsigned int j = 0; j < this_block_size; j++)
1176 irregular[counter_macro++] =
1177 irregular_cells[present_block *
block_size + j];
1182 irregular_cells.swap(irregular);
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)
1207 const std::vector<unsigned int> &cell_active_fe_index,
1209 std::vector<unsigned int> & renumbering,
1210 std::vector<unsigned char> & irregular_cells,
1214 if (n_macro_cells == 0)
1224 connectivity_blocks);
1227 if (
scheme == partition_color ||
1236 std::vector<unsigned int> cell_partition(
n_blocks,
1242 std::vector<unsigned int> partition_list(
n_blocks, 0);
1243 std::vector<unsigned int> partition_2layers_list(
n_blocks, 0);
1246 std::vector<unsigned int> partition_size(2, 0);
1248 unsigned int partition = 0;
1254 unsigned int cluster_size = 1;
1255 if (
scheme == partition_partition)
1275 if (
scheme == partition_partition)
1281 cell_active_fe_index,
1288 partition_2layers_list,
1298 partition_2layers_list);
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)
1313 renumbering_in.swap(renumbering);
1314 if (
scheme == partition_partition)
1319 for (
unsigned int j = 0; j < renumbering.size(); j++)
1320 renumbering[j] = renumbering_in[partition_2layers_list[j]];
1329 std::vector<unsigned int> block_start(n_macro_cells + 1);
1330 std::vector<unsigned char> irregular(n_macro_cells);
1332 unsigned int counter = 0;
1333 unsigned int mcell_start = 0;
1335 for (
unsigned int block = 0; block <
n_blocks; block++)
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);
1342 unsigned int n_comp = (irregular_cells[mcell] > 0) ?
1343 irregular_cells[mcell] :
1345 block_start[block + 1] += n_comp;
1351 unsigned int counter_macro = 0;
1352 unsigned int block_size_last =
1354 if (block_size_last == 0)
1357 unsigned int tick = 0;
1358 for (
unsigned int block = 0; block <
n_blocks; block++)
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];
1364 renumbering[counter++] = renumbering_in[cell];
1365 unsigned int this_block_size =
1373 for (
unsigned int j = 0; j < this_block_size; j++)
1374 irregular[counter_macro++] =
1375 irregular_cells[present_block *
block_size + j];
1380 irregular_cells.swap(irregular);
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)
1402 const std::vector<unsigned int> &cell_active_fe_index,
1404 std::vector<unsigned int> & renumbering,
1405 std::vector<unsigned char> & irregular_cells,
1409 if (n_macro_cells == 0)
1426 std::vector<unsigned int> partition_partition_list(
n_active_cells, 0);
1429 std::vector<unsigned int> partition_size(2, 0);
1431 unsigned int partition = 0;
1445 cell_active_fe_index,
1452 partition_partition_list,
1455 partition_list.swap(renumbering);
1457 for (
unsigned int j = 0; j < renumbering.size(); j++)
1458 renumbering[j] = partition_list[partition_partition_list[j]];
1470 const std::vector<unsigned char> &irregular_cells,
1474 std::vector<std::vector<unsigned int>> cell_blocks(
n_blocks);
1476 unsigned int cell = 0;
1477 for (
unsigned int i = 0, mcell = 0; i <
n_blocks; ++i)
1479 for (
unsigned int c = 0;
1483 unsigned int ncomp = (irregular_cells[mcell] > 0) ?
1484 irregular_cells[mcell] :
1486 for (
unsigned int c = 0; c < ncomp; ++c, ++cell)
1488 cell_blocks[i].push_back(cell);
1489 touched_cells[cell] = i;
1494 for (
unsigned int i = 0; i < cell_blocks.size(); ++i)
1495 for (
unsigned int col = 0; col < cell_blocks[i].size(); ++col)
1498 connectivity_cells.
begin(cell_blocks[i][col]);
1499 it != connectivity_cells.
end(cell_blocks[i][col]);
1502 if (touched_cells[it->column()] != i)
1503 connectivity_blocks.
add(i, touched_cells[it->column()]);
1515 const std::vector<unsigned int> &cell_active_fe_index,
1516 const unsigned int partition,
1517 const unsigned int cluster_size,
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)
1526 const unsigned int n_ghost_slots =
1530 std::vector<unsigned int> neighbor_list;
1533 std::vector<unsigned int> neighbor_neighbor_list;
1537 irregular_cells.back() = 0;
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);
1548 unsigned int n_macro_cells_before = 0;
1554 std::vector<unsigned int> cell_partition_l2(
1560 unsigned int counter = 0;
1561 unsigned int missing_macros;
1562 for (
unsigned int part = 0; part < partition; ++part)
1564 neighbor_neighbor_list.resize(0);
1565 neighbor_list.resize(0);
1567 unsigned int partition_l2 = 0;
1568 unsigned int start_up = partition_size[part];
1569 unsigned int partition_counter = 0;
1572 if (neighbor_list.size() == 0)
1575 partition_counter = 0;
1576 for (
unsigned int j = start_up;
1577 j < partition_size[part + 1];
1579 if (cell_partition[partition_list[j]] == part &&
1580 cell_partition_l2[partition_list[j]] ==
1585 partition_counter = 1;
1589 cell_partition_l2[partition_list[start_up]] =
1591 neighbor_neighbor_list.push_back(
1592 partition_list[start_up]);
1593 partition_partition_list[counter++] =
1594 partition_list[start_up];
1601 partition_counter = 0;
1602 for (
const unsigned int neighbor : neighbor_list)
1604 Assert(cell_partition[neighbor] == part,
1606 Assert(cell_partition_l2[neighbor] == partition_l2 - 1,
1608 auto neighbor_it = connectivity.
begin(neighbor);
1609 const auto end_it = connectivity.
end(neighbor);
1610 for (; neighbor_it != end_it; ++neighbor_it)
1612 if (cell_partition[neighbor_it->column()] == part &&
1613 cell_partition_l2[neighbor_it->column()] ==
1616 cell_partition_l2[neighbor_it->column()] =
1618 neighbor_neighbor_list.push_back(
1619 neighbor_it->column());
1620 partition_partition_list[counter++] =
1621 neighbor_it->column();
1622 partition_counter++;
1627 if (partition_counter > 0)
1629 int index_before = neighbor_neighbor_list.size(),
1630 index = index_before;
1635 std::vector<unsigned int> remaining_per_macro_cell(
1637 std::vector<std::vector<unsigned int>>
1638 renumbering_fe_index;
1641 if (hp_bool ==
true)
1643 renumbering_fe_index.resize(max_fe_index + 1);
1644 for (cell = counter - partition_counter;
1648 renumbering_fe_index
1649 [cell_active_fe_index.empty() ?
1651 cell_active_fe_index
1652 [partition_partition_list[cell]]]
1653 .push_back(partition_partition_list[cell]);
1656 for (
unsigned int j = 0; j < max_fe_index + 1; j++)
1658 remaining_per_macro_cell[j] =
1659 renumbering_fe_index[j].size() %
1661 if (remaining_per_macro_cell[j] != 0)
1664 ((renumbering_fe_index[j].size() +
1671 remaining_per_macro_cell.resize(1);
1672 remaining_per_macro_cell[0] =
1676 if (remaining_per_macro_cell[0] != 0)
1683 cluster_size - (missing_macros % cluster_size);
1686 while (missing_macros > 0 || filled ==
false)
1690 index = neighbor_neighbor_list.size();
1691 if (index == index_before)
1693 if (missing_macros != 0)
1695 neighbor_neighbor_list.resize(0);
1700 index_before = index;
1703 unsigned int additional =
1704 neighbor_neighbor_list[index];
1715 for (; neighbor != end; ++neighbor)
1717 if (cell_partition[neighbor->
column()] == part &&
1718 cell_partition_l2[neighbor->
column()] ==
1721 unsigned int this_index = 0;
1722 if (hp_bool ==
true)
1724 cell_active_fe_index.empty() ?
1726 cell_active_fe_index[neighbor
1733 if (missing_macros > 0 ||
1734 remaining_per_macro_cell[this_index] > 0)
1736 cell_partition_l2[neighbor->
column()] =
1738 neighbor_neighbor_list.push_back(
1740 if (hp_bool ==
true)
1741 renumbering_fe_index[this_index]
1742 .push_back(neighbor->
column());
1743 partition_partition_list[counter] =
1746 partition_counter++;
1747 if (remaining_per_macro_cell
1748 [this_index] == 0 &&
1751 remaining_per_macro_cell[this_index]++;
1752 if (remaining_per_macro_cell
1756 remaining_per_macro_cell[this_index] =
1759 if (missing_macros == 0)
1762 for (
unsigned int fe_ind = 0;
1763 fe_ind < max_fe_index + 1;
1765 if (remaining_per_macro_cell
1775 if (hp_bool ==
true)
1780 cell = counter - partition_counter;
1781 for (
unsigned int j = 0; j < max_fe_index + 1; j++)
1783 for (
const unsigned int jj :
1784 renumbering_fe_index[j])
1785 renumbering[cell++] = jj;
1786 if (renumbering_fe_index[j].size() %
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() +
1798 renumbering_fe_index[j].resize(0);
1803 n_macro_cells_before +=
1807 irregular_cells[n_macro_cells_before] =
1809 n_macro_cells_before++;
1816 neighbor_list = neighbor_neighbor_list;
1817 neighbor_neighbor_list.resize(0);
1823 if (hp_bool ==
true)
1825 partition_partition_list.swap(renumbering);
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)
1843 std::vector<unsigned int> cell_color(
n_blocks, n_macro_cells);
1844 std::vector<bool> color_finder;
1848 unsigned int color_counter = 0, index_counter = 0;
1849 for (
unsigned int part = 0; part < partition; part++)
1852 unsigned int max_color = 0;
1853 for (
unsigned int k = partition_size[part];
1854 k < partition_size[part + 1];
1857 unsigned int cell = partition_list[k];
1858 unsigned int n_neighbors = connectivity.
row_length(cell);
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)
1872 if (cell_partition[neighbor->
column()] == part &&
1873 cell_color[neighbor->
column()] <= n_neighbors)
1874 color_finder[cell_color[neighbor->
column()]] =
false;
1877 cell_color[cell] = 0;
1878 while (color_finder[cell_color[cell]] ==
false)
1880 if (cell_color[cell] > max_color)
1881 max_color = cell_color[cell];
1886 for (
unsigned int color = 0; color <= max_color; color++)
1890 for (
unsigned int k = partition_size[part];
1891 k < partition_size[part + 1];
1894 unsigned int cell = partition_list[k];
1895 if (cell_color[cell] == color)
1897 partition_color_list[color_counter++] = cell;
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 1924 std::vector<unsigned int> neighbor_list;
1927 std::vector<unsigned int> neighbor_neighbor_list;
1937 unsigned int counter = 0;
1938 unsigned int start_nonboundary =
1945 if (n_macro_cells == 0)
1948 start_nonboundary = n_macro_cells;
1949 if (
scheme == partition_color ||
1953 if (
scheme == partition_color ||
1963 unsigned int start_up = 0;
1965 unsigned int remainder = cluster_size;
1973 if (start_nonboundary > 0)
1975 for (
unsigned int cell = 0; cell < start_nonboundary; ++cell)
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()++;
1983 start_nonboundary = 0;
1984 remainder -= (start_nonboundary % cluster_size);
1985 if (remainder == cluster_size)
1992 cell_partition[start_up] = partition;
1993 neighbor_list.push_back(start_up);
1994 partition_list[counter++] = start_up;
1995 partition_size.back()++;
1998 if (remainder == cluster_size)
2001 int index_before = neighbor_list.size(), index = index_before,
2003 while (remainder > 0)
2005 if (index == index_stop)
2007 index = neighbor_list.size();
2008 if (index == index_before)
2010 neighbor_list.resize(0);
2013 index_stop = index_before;
2014 index_before = index;
2017 unsigned int additional = neighbor_list[index];
2019 connectivity.
begin(additional),
2021 connectivity.
end(additional);
2022 for (; neighbor != end; ++neighbor)
2024 if (cell_partition[neighbor->
column()] ==
2027 partition_size.back()++;
2028 cell_partition[neighbor->
column()] = partition;
2029 neighbor_list.push_back(neighbor->
column());
2030 partition_list[counter++] = neighbor->
column();
2038 while (neighbor_list.size() > 0)
2043 unsigned int partition_counter = 0;
2046 partition_size.push_back(partition_size.back());
2050 for (
const unsigned int cell : neighbor_list)
2052 Assert(cell_partition[cell] == partition - 1,
2054 auto neighbor = connectivity.
begin(cell);
2055 const auto end = connectivity.
end(cell);
2056 for (; neighbor != end; ++neighbor)
2058 if (cell_partition[neighbor->column()] ==
2061 partition_size.back()++;
2062 cell_partition[neighbor->column()] = partition;
2066 neighbor_neighbor_list.push_back(neighbor->column());
2067 partition_list[counter++] = neighbor->column();
2068 partition_counter++;
2072 remainder = cluster_size - (partition_counter % cluster_size);
2073 if (remainder == cluster_size)
2076 int index_before = neighbor_neighbor_list.size(),
2077 index = index_before;
2078 while (remainder > 0)
2080 if (index == index_stop)
2082 index = neighbor_neighbor_list.size();
2083 if (index == index_before)
2085 neighbor_neighbor_list.resize(0);
2088 index_stop = index_before;
2089 index_before = index;
2092 unsigned int additional = neighbor_neighbor_list[index];
2096 end = connectivity.
end(
2098 for (; neighbor != end; ++neighbor)
2100 if (cell_partition[neighbor->
column()] ==
2103 partition_size.back()++;
2104 cell_partition[neighbor->
column()] = partition;
2105 neighbor_neighbor_list.push_back(neighbor->
column());
2106 partition_list[counter++] = neighbor->
column();
2114 neighbor_list = neighbor_neighbor_list;
2115 neighbor_neighbor_list.resize(0);
2121 for (
unsigned int j = start_up; j <
n_blocks; ++j)
2127 remainder = cluster_size;
2141 evens = (partition + 1) / 2;
2142 odds = partition / 2;
2150 for (
unsigned int part = 0; part < partition; part++)
2171 internal::MatrixFreeFunctions::TaskInfo::print_memory_statistics<std::ostream>(
2173 const std::size_t)
const;
2179 DEAL_II_NAMESPACE_CLOSE
static const unsigned int invalid_unsigned_int
#define AssertDimension(dim1, dim2)
virtual void cell(const std::pair< unsigned int, unsigned int > &cell_range)=0
TasksParallelScheme scheme
void parallel_for(Iterator x_begin, Iterator x_end, const Functor &functor, const unsigned int grainsize)
void loop(MFWorkerInterface &worker) const
void add(const size_type i, const size_type j)
void guess_block_size(const unsigned int dofs_per_cell)
std::vector< unsigned int > partition_row_index
#define AssertIndexRange(index, range)
std::vector< unsigned int > partition_n_workers
virtual void vector_update_ghosts_finish()=0
Finishes the communication for the update ghost values operation.
std::vector< unsigned int > cell_partition_data
#define AssertThrow(cond, exc)
void make_connectivity_cells_to_blocks(const std::vector< unsigned char > &irregular_cells, const DynamicSparsityPattern &connectivity_cells, DynamicSparsityPattern &connectivity_blocks) const
std::vector< unsigned int > partition_n_blocked_workers
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)
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)
unsigned int vectorization_length
#define Assert(cond, exc)
unsigned int n_blocked_workers
void make_thread_graph_partition_color(DynamicSparsityPattern &connectivity, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &irregular_cells, const bool hp_bool)
void update_task_info(const unsigned int partition)
std::vector< unsigned int > boundary_partition_data
void initial_setup_blocks_tasks(const std::vector< unsigned int > &boundary_cells, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &incompletely_filled_vectorization)
unsigned int n_ghost_cells
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)
std::vector< unsigned int > invert_permutation(const std::vector< unsigned int > &permutation)
unsigned int n_active_cells
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
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)
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)
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
virtual void vector_compress_start()=0
Starts the communication for the vector compress operation.
std::vector< unsigned int > partition_odds
std::size_t memory_consumption() const
static ::ExceptionBase & ExcNotImplemented()
static unsigned int n_threads()
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
std::vector< unsigned int > partition_evens
void print_memory_statistics(StreamType &out, std::size_t data_length) const
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)
virtual void vector_compress_finish()=0
Finishes the communication for the vector compress operation.