45 #define BOOST_BIND_GLOBAL_PLACEHOLDERS
46 #include <boost/config.hpp>
47 #include <boost/graph/adjacency_list.hpp>
48 #include <boost/graph/bandwidth.hpp>
49 #include <boost/graph/cuthill_mckee_ordering.hpp>
50 #include <boost/graph/king_ordering.hpp>
51 #include <boost/graph/minimum_degree_ordering.hpp>
52 #include <boost/graph/properties.hpp>
53 #include <boost/random.hpp>
54 #include <boost/random/uniform_int_distribution.hpp>
55 #undef BOOST_BIND_GLOBAL_PLACEHOLDERS
73 using namespace ::
boost;
76 using Graph = adjacency_list<vecS,
79 property<vertex_color_t,
81 property<vertex_degree_t, int>>>;
82 using Vertex = graph_traits<Graph>::vertex_descriptor;
83 using size_type = graph_traits<Graph>::vertices_size_type;
85 using Pair = std::pair<size_type, size_type>;
91 template <
typename DoFHandlerType>
94 const bool use_constraints,
97 boosttypes::vertex_degree_t>::type
108 dof_handler.n_dofs());
112 for (
unsigned int row = 0; row < dsp.
n_rows(); ++row)
113 for (
unsigned int col = 0; col < dsp.
row_length(row); ++col)
117 boosttypes::graph_traits<boosttypes::Graph>::vertex_iterator ui, ui_end;
119 graph_degree = get(::boost::vertex_degree, graph);
120 for (::boost::tie(ui, ui_end) =
vertices(graph); ui != ui_end; ++ui)
121 graph_degree[*ui] = degree(*ui, graph);
126 template <
typename DoFHandlerType>
129 const bool reversed_numbering,
130 const bool use_constraints)
132 std::vector<types::global_dof_index> renumbering(
142 dof_handler.renumber_dofs(renumbering);
146 template <
typename DoFHandlerType>
149 const DoFHandlerType & dof_handler,
150 const bool reversed_numbering,
151 const bool use_constraints)
155 boosttypes::vertex_degree_t>::type graph_degree;
160 boosttypes::vertex_index_t>::type index_map =
161 get(::boost::vertex_index, graph);
164 std::vector<boosttypes::Vertex> inv_perm(num_vertices(graph));
166 if (reversed_numbering ==
false)
167 ::boost::cuthill_mckee_ordering(graph,
169 get(::boost::vertex_color, graph),
170 make_degree_map(graph));
172 ::boost::cuthill_mckee_ordering(graph,
174 get(::boost::vertex_color, graph),
175 make_degree_map(graph));
178 new_dof_indices[index_map[inv_perm[c]]] = c;
180 Assert(std::find(new_dof_indices.begin(),
181 new_dof_indices.end(),
188 template <
typename DoFHandlerType>
191 const bool reversed_numbering,
192 const bool use_constraints)
194 std::vector<types::global_dof_index> renumbering(
204 dof_handler.renumber_dofs(renumbering);
208 template <
typename DoFHandlerType>
211 const DoFHandlerType & dof_handler,
212 const bool reversed_numbering,
213 const bool use_constraints)
217 boosttypes::vertex_degree_t>::type graph_degree;
222 boosttypes::vertex_index_t>::type index_map =
223 get(::boost::vertex_index, graph);
226 std::vector<boosttypes::Vertex> inv_perm(num_vertices(graph));
228 if (reversed_numbering ==
false)
234 new_dof_indices[index_map[inv_perm[c]]] = c;
236 Assert(std::find(new_dof_indices.begin(),
237 new_dof_indices.end(),
244 template <
typename DoFHandlerType>
247 const bool reversed_numbering,
248 const bool use_constraints)
250 std::vector<types::global_dof_index> renumbering(
260 dof_handler.renumber_dofs(renumbering);
264 template <
typename DoFHandlerType>
267 std::vector<types::global_dof_index> &new_dof_indices,
268 const DoFHandlerType & dof_handler,
269 const bool reversed_numbering,
270 const bool use_constraints)
272 (void)use_constraints;
281 using namespace ::
boost;
286 using Graph = adjacency_list<vecS, vecS, directedS>;
288 int n = dof_handler.n_dofs();
292 std::vector<::types::global_dof_index> dofs_on_this_cell;
294 typename DoFHandlerType::active_cell_iterator cell = dof_handler
296 endc = dof_handler.end();
298 for (; cell != endc; ++cell)
300 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
302 dofs_on_this_cell.resize(dofs_per_cell);
304 cell->get_active_or_mg_dof_indices(dofs_on_this_cell);
305 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
306 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
307 if (dofs_on_this_cell[i] > dofs_on_this_cell[j])
309 add_edge(dofs_on_this_cell[i], dofs_on_this_cell[j], G);
310 add_edge(dofs_on_this_cell[j], dofs_on_this_cell[i], G);
315 using Vector = std::vector<int>;
318 Vector inverse_perm(n, 0);
323 Vector supernode_sizes(n, 1);
326 ::boost::property_map<Graph, vertex_index_t>::type
id =
327 get(vertex_index, G);
333 minimum_degree_ordering(
335 make_iterator_property_map(degree.begin(),
id, degree[0]),
338 make_iterator_property_map(supernode_sizes.begin(),
345 for (
int i = 0; i < n; ++i)
347 Assert(std::find(perm.begin(), perm.end(), i) != perm.end(),
349 Assert(std::find(inverse_perm.begin(), inverse_perm.end(), i) !=
355 if (reversed_numbering ==
true)
356 std::copy(perm.begin(), perm.end(), new_dof_indices.begin());
360 new_dof_indices.begin());
367 template <
typename DoFHandlerType>
370 const bool reversed_numbering,
371 const bool use_constraints,
372 const std::vector<types::global_dof_index> &starting_indices)
374 std::vector<types::global_dof_index> renumbering(
375 dof_handler.locally_owned_dofs().n_elements(),
386 dof_handler.renumber_dofs(renumbering);
391 template <
typename DoFHandlerType>
394 std::vector<types::global_dof_index> & new_indices,
395 const DoFHandlerType & dof_handler,
396 const bool reversed_numbering,
397 const bool use_constraints,
398 const std::vector<types::global_dof_index> &starting_indices)
402 if (dof_handler.locally_owned_dofs().n_elements() == 0)
418 constraints.
reinit(locally_relevant_dofs);
423 const IndexSet &locally_owned_dofs = dof_handler.locally_owned_dofs();
436 if (reversed_numbering)
452 bool needs_locally_active =
false;
453 for (
const auto starting_index : starting_indices)
455 if ((needs_locally_active ==
457 (locally_owned_dofs.
is_element(starting_index) ==
false))
460 locally_active_dofs.
is_element(starting_index),
462 "You specified global degree of freedom " +
464 " as a starting index, but this index is not among the "
465 "locally active ones on this processor, as required "
466 "for this function."));
467 needs_locally_active =
true;
472 (needs_locally_active ? locally_active_dofs : locally_owned_dofs);
478 dof_handler.n_dofs(),
484 std::vector<types::global_dof_index> row_entries;
485 for (
unsigned int i = 0; i < index_set_to_use.
n_elements(); ++i)
489 const unsigned int row_length = dsp.
row_length(row);
491 for (
unsigned int j = 0; j < row_length; ++j)
494 if (col != row && index_set_to_use.
is_element(col))
504 std::vector<types::global_dof_index> local_starting_indices(
505 starting_indices.size());
506 for (
unsigned int i = 0; i < starting_indices.size(); ++i)
507 local_starting_indices[i] =
512 std::vector<types::global_dof_index> my_new_indices(
516 local_starting_indices);
517 if (reversed_numbering)
526 if (needs_locally_active ==
true)
529 IndexSet active_but_not_owned_dofs = locally_active_dofs;
530 active_but_not_owned_dofs.
subtract_set(locally_owned_dofs);
532 std::set<types::global_dof_index> erase_these_indices;
533 for (
const auto p : active_but_not_owned_dofs)
538 erase_these_indices.insert(my_new_indices[index]);
541 Assert(erase_these_indices.size() ==
542 active_but_not_owned_dofs.n_elements(),
544 Assert(
static_cast<unsigned int>(
545 std::count(my_new_indices.begin(),
546 my_new_indices.end(),
548 active_but_not_owned_dofs.n_elements(),
552 std::vector<types::global_dof_index> translate_indices(
553 my_new_indices.size());
555 std::set<types::global_dof_index>::const_iterator
556 next_erased_index = erase_these_indices.begin();
558 for (
unsigned int i = 0; i < translate_indices.size(); ++i)
559 if ((next_erased_index != erase_these_indices.end()) &&
560 (*next_erased_index == i))
567 translate_indices[i] = next_new_index;
577 new_indices.reserve(locally_owned_dofs.
n_elements());
578 for (
const auto &p : my_new_indices)
583 new_indices.push_back(translate_indices[p]);
589 new_indices = std::move(my_new_indices);
602 template <
typename DoFHandlerType>
605 const unsigned int level,
606 const bool reversed_numbering,
607 const std::vector<types::global_dof_index> &starting_indices)
614 dof_handler.n_dofs(
level));
617 std::vector<types::global_dof_index> new_indices(dsp.
n_rows());
620 if (reversed_numbering)
626 dof_handler.renumber_dofs(
level, new_indices);
631 template <
typename DoFHandlerType>
634 const std::vector<unsigned int> &component_order_arg)
636 std::vector<types::global_dof_index> renumbering(
639 typename DoFHandlerType::active_cell_iterator start =
640 dof_handler.begin_active();
641 const typename DoFHandlerType::level_cell_iterator
end = dof_handler.end();
645 DoFHandlerType::space_dimension>(
646 renumbering, start,
end, component_order_arg,
false);
657 Assert((result == dof_handler.n_locally_owned_dofs()) ||
658 ((dof_handler.n_locally_owned_dofs() < dof_handler.n_dofs()) &&
659 (result <= dof_handler.n_dofs())),
662 dof_handler.renumber_dofs(renumbering);
667 template <
typename DoFHandlerType>
670 const unsigned int level,
671 const std::vector<unsigned int> &component_order_arg)
676 std::vector<types::global_dof_index> renumbering(
677 dof_handler.locally_owned_mg_dofs(
level).n_elements(),
680 typename DoFHandlerType::level_cell_iterator start =
681 dof_handler.begin(
level);
682 typename DoFHandlerType::level_cell_iterator
end = dof_handler.end(
level);
686 DoFHandlerType::space_dimension>(
687 renumbering, start,
end, component_order_arg,
true);
694 if (renumbering.size() != 0)
695 dof_handler.renumber_dofs(
level, renumbering);
700 template <
int dim,
int spacedim,
typename CellIterator>
703 const CellIterator & start,
705 const std::vector<unsigned int> &component_order_arg,
706 const bool is_level_operation)
709 start->get_dof_handler().get_fe_collection());
715 new_indices.resize(0);
722 const IndexSet &locally_owned_dofs =
724 start->get_dof_handler().locally_owned_mg_dofs(start->level()) :
725 start->get_dof_handler().locally_owned_dofs();
729 std::vector<unsigned int> component_order(component_order_arg);
734 if (component_order.size() == 0)
735 for (
unsigned int i = 0; i < fe_collection.
n_components(); ++i)
736 component_order.push_back(i);
742 for (
const unsigned int component : component_order)
750 std::vector<types::global_dof_index> local_dof_indices;
762 std::vector<std::vector<unsigned int>> component_list(fe_collection.
size());
763 for (
unsigned int f = 0; f < fe_collection.
size(); ++f)
767 component_list[f].resize(dofs_per_cell);
768 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
770 component_list[f][i] =
774 const unsigned int comp =
780 component_list[f][i] = component_order[comp];
795 std::vector<std::vector<types::global_dof_index>> component_to_dof_map(
797 for (CellIterator cell = start; cell !=
end; ++cell)
799 if (is_level_operation)
802 if (!cell->is_locally_owned_on_level())
809 if (!cell->is_locally_owned())
813 if (is_level_operation)
815 cell->level() == start->level(),
817 "Multigrid renumbering in compute_component_wise() needs to be applied to a single level!"));
822 const unsigned int fe_index = cell->active_fe_index();
823 const unsigned int dofs_per_cell =
824 fe_collection[fe_index].dofs_per_cell;
825 local_dof_indices.resize(dofs_per_cell);
826 cell->get_active_or_mg_dof_indices(local_dof_indices);
828 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
829 if (locally_owned_dofs.
is_element(local_dof_indices[i]))
830 component_to_dof_map[component_list[fe_index][i]].push_back(
831 local_dof_indices[i]);
857 for (
unsigned int component = 0; component < fe_collection.
n_components();
860 std::sort(component_to_dof_map[component].
begin(),
861 component_to_dof_map[component].
end());
862 component_to_dof_map[component].erase(
863 std::unique(component_to_dof_map[component].
begin(),
864 component_to_dof_map[component].
end()),
865 component_to_dof_map[component].
end());
870 const unsigned int n_buckets = fe_collection.
n_components();
871 std::vector<types::global_dof_index> shifts(n_buckets);
875 &start->get_dof_handler().get_triangulation())))
877 #ifdef DEAL_II_WITH_MPI
878 std::vector<types::global_dof_index> local_dof_count(n_buckets);
880 for (
unsigned int c = 0; c < n_buckets; ++c)
881 local_dof_count[c] = component_to_dof_map[c].size();
883 std::vector<types::global_dof_index> prefix_dof_count(n_buckets);
884 const int ierr = MPI_Exscan(local_dof_count.data(),
885 prefix_dof_count.data(),
889 tria->get_communicator());
892 std::vector<types::global_dof_index> global_dof_count(n_buckets);
894 tria->get_communicator(),
899 for (
unsigned int c = 0; c < n_buckets; ++c)
901 shifts[c] = prefix_dof_count[c] + cumulated;
902 cumulated += global_dof_count[c];
912 for (
unsigned int c = 1; c < fe_collection.
n_components(); ++c)
913 shifts[c] = shifts[c - 1] + component_to_dof_map[c - 1].size();
922 for (
unsigned int component = 0; component < fe_collection.
n_components();
925 const typename std::vector<types::global_dof_index>::const_iterator
926 begin_of_component = component_to_dof_map[component].begin(),
927 end_of_component = component_to_dof_map[component].end();
929 next_free_index = shifts[component];
931 for (
typename std::vector<types::global_dof_index>::const_iterator
932 dof_index = begin_of_component;
933 dof_index != end_of_component;
944 return next_free_index;
949 template <
int dim,
int spacedim>
953 std::vector<types::global_dof_index> renumbering(
981 (result <= dof_handler.
n_dofs())),
989 template <
int dim,
int spacedim>
993 std::vector<types::global_dof_index> renumbering(
1020 template <
int dim,
int spacedim>
1027 std::vector<types::global_dof_index> renumbering(
1050 if (renumbering.size() != 0)
1056 template <
int dim,
int spacedim,
class ITERATOR,
class ENDITERATOR>
1059 const ITERATOR & start,
1060 const ENDITERATOR &
end,
1061 const bool is_level_operation)
1064 start->get_dof_handler().get_fe_collection());
1070 new_indices.resize(0);
1077 const IndexSet &locally_owned_dofs =
1078 is_level_operation ?
1079 start->get_dof_handler().locally_owned_mg_dofs(start->level()) :
1080 start->get_dof_handler().locally_owned_dofs();
1084 std::vector<types::global_dof_index> local_dof_indices;
1089 std::vector<std::vector<types::global_dof_index>> block_list(
1090 fe_collection.
size());
1091 for (
unsigned int f = 0; f < fe_collection.
size(); ++f)
1110 std::vector<std::vector<types::global_dof_index>> block_to_dof_map(
1112 for (ITERATOR cell = start; cell !=
end; ++cell)
1114 if (is_level_operation)
1117 if (!cell->is_locally_owned_on_level())
1124 if (!cell->is_locally_owned())
1128 if (is_level_operation)
1130 cell->level() == start->level(),
1132 "Multigrid renumbering in compute_block_wise() needs to be applied to a single level!"));
1137 const unsigned int fe_index = cell->active_fe_index();
1138 const unsigned int dofs_per_cell =
1139 fe_collection[fe_index].dofs_per_cell;
1140 local_dof_indices.resize(dofs_per_cell);
1141 cell->get_active_or_mg_dof_indices(local_dof_indices);
1143 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1144 if (locally_owned_dofs.
is_element(local_dof_indices[i]))
1145 block_to_dof_map[block_list[fe_index][i]].push_back(
1146 local_dof_indices[i]);
1161 for (
unsigned int block = 0; block < fe_collection.
n_blocks(); ++block)
1163 std::sort(block_to_dof_map[block].
begin(),
1164 block_to_dof_map[block].
end());
1165 block_to_dof_map[block].erase(
1166 std::unique(block_to_dof_map[block].
begin(),
1167 block_to_dof_map[block].
end()),
1168 block_to_dof_map[block].
end());
1173 const unsigned int n_buckets = fe_collection.
n_blocks();
1174 std::vector<types::global_dof_index> shifts(n_buckets);
1178 &start->get_dof_handler().get_triangulation())))
1180 #ifdef DEAL_II_WITH_MPI
1181 std::vector<types::global_dof_index> local_dof_count(n_buckets);
1183 for (
unsigned int c = 0; c < n_buckets; ++c)
1184 local_dof_count[c] = block_to_dof_map[c].size();
1186 std::vector<types::global_dof_index> prefix_dof_count(n_buckets);
1187 const int ierr = MPI_Exscan(local_dof_count.data(),
1188 prefix_dof_count.data(),
1192 tria->get_communicator());
1195 std::vector<types::global_dof_index> global_dof_count(n_buckets);
1197 tria->get_communicator(),
1202 for (
unsigned int c = 0; c < n_buckets; ++c)
1204 shifts[c] = prefix_dof_count[c] + cumulated;
1205 cumulated += global_dof_count[c];
1215 for (
unsigned int c = 1; c < fe_collection.
n_blocks(); ++c)
1216 shifts[c] = shifts[c - 1] + block_to_dof_map[c - 1].size();
1225 for (
unsigned int block = 0; block < fe_collection.
n_blocks(); ++block)
1227 const typename std::vector<types::global_dof_index>::const_iterator
1228 begin_of_component = block_to_dof_map[block].begin(),
1229 end_of_component = block_to_dof_map[block].end();
1231 next_free_index = shifts[block];
1233 for (
typename std::vector<types::global_dof_index>::const_iterator
1234 dof_index = begin_of_component;
1235 dof_index != end_of_component;
1246 return next_free_index;
1258 template <
int dim,
class CellIteratorType>
1260 compute_hierarchical_recursive(
1263 const CellIteratorType & cell,
1264 const IndexSet & locally_owned_dof_indices,
1265 std::vector<types::global_dof_index> &new_indices)
1268 next_free_dof_offset;
1270 if (cell->has_children())
1273 for (
unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell;
1275 current_next_free_dof_offset =
1276 compute_hierarchical_recursive<dim>(current_next_free_dof_offset,
1279 locally_owned_dof_indices,
1299 if (cell->is_locally_owned())
1302 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
1303 std::vector<types::global_dof_index> local_dof_indices(
1305 cell->get_dof_indices(local_dof_indices);
1324 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1325 if (locally_owned_dof_indices.
is_element(local_dof_indices[i]))
1329 const unsigned int idx =
1331 local_dof_indices[i]);
1335 my_starting_index + current_next_free_dof_offset;
1336 ++current_next_free_dof_offset;
1342 return current_next_free_dof_offset;
1348 template <
typename DoFHandlerType>
1352 const int dim = DoFHandlerType::dimension;
1353 const int spacedim = DoFHandlerType::space_dimension;
1355 std::vector<types::global_dof_index> renumbering(
1359 const IndexSet locally_owned = dof_handler.locally_owned_dofs();
1376 &dof_handler.get_triangulation()))
1378 #ifdef DEAL_II_WITH_MPI
1380 dof_handler.locally_owned_dofs().n_elements();
1381 MPI_Exscan(&local_size,
1386 tria->get_communicator());
1392 *
>(&dof_handler.get_triangulation()))
1394 #ifdef DEAL_II_WITH_P4EST
1399 for (
unsigned int c = 0; c < tria->n_cells(0); ++c)
1401 const unsigned int coarse_cell_index =
1402 tria->get_p4est_tree_to_coarse_cell_permutation()[c];
1404 const typename DoFHandlerType::level_cell_iterator this_cell(
1405 tria, 0, coarse_cell_index, &dof_handler);
1407 next_free_dof_offset =
1408 compute_hierarchical_recursive<dim>(next_free_dof_offset,
1422 for (
typename DoFHandlerType::cell_iterator cell = dof_handler.begin(0);
1423 cell != dof_handler.end(0);
1425 next_free_dof_offset =
1426 compute_hierarchical_recursive<dim>(next_free_dof_offset,
1437 Assert((next_free_dof_offset == dof_handler.n_locally_owned_dofs()) ||
1438 ((dof_handler.n_locally_owned_dofs() < dof_handler.n_dofs()) &&
1439 (next_free_dof_offset <= dof_handler.n_dofs())),
1443 Assert(std::find(renumbering.begin(),
1448 dof_handler.renumber_dofs(renumbering);
1453 template <
typename DoFHandlerType>
1456 const std::vector<bool> &selected_dofs)
1458 std::vector<types::global_dof_index> renumbering(
1462 dof_handler.renumber_dofs(renumbering);
1467 template <
typename DoFHandlerType>
1470 const std::vector<bool> &selected_dofs,
1471 const unsigned int level)
1476 std::vector<types::global_dof_index> renumbering(
1483 dof_handler.renumber_dofs(
level, renumbering);
1488 template <
typename DoFHandlerType>
1491 std::vector<types::global_dof_index> &new_indices,
1492 const DoFHandlerType & dof_handler,
1493 const std::vector<bool> & selected_dofs)
1496 Assert(selected_dofs.size() == n_dofs,
1501 Assert(new_indices.size() == n_dofs,
1505 std::count(selected_dofs.begin(), selected_dofs.end(),
false);
1510 if (selected_dofs[i] ==
false)
1512 new_indices[i] = next_unselected;
1517 new_indices[i] = next_selected;
1526 template <
typename DoFHandlerType>
1529 std::vector<types::global_dof_index> &new_indices,
1530 const DoFHandlerType & dof_handler,
1531 const std::vector<bool> & selected_dofs,
1532 const unsigned int level)
1537 const unsigned int n_dofs = dof_handler.n_dofs(
level);
1538 Assert(selected_dofs.size() == n_dofs,
1543 Assert(new_indices.size() == n_dofs,
1546 const unsigned int n_selected_dofs =
1547 std::count(selected_dofs.begin(), selected_dofs.end(),
false);
1549 unsigned int next_unselected = 0;
1550 unsigned int next_selected = n_selected_dofs;
1551 for (
unsigned int i = 0; i < n_dofs; ++i)
1552 if (selected_dofs[i] ==
false)
1554 new_indices[i] = next_unselected;
1559 new_indices[i] = next_selected;
1568 template <
typename DoFHandlerType>
1571 DoFHandlerType & dof,
1572 const std::vector<typename DoFHandlerType::active_cell_iterator> &cells)
1574 std::vector<types::global_dof_index> renumbering(
1575 dof.n_locally_owned_dofs());
1576 std::vector<types::global_dof_index> reverse(dof.n_locally_owned_dofs());
1579 dof.renumber_dofs(renumbering);
1583 template <
typename DoFHandlerType>
1586 std::vector<types::global_dof_index> &new_indices,
1587 std::vector<types::global_dof_index> &reverse,
1588 const DoFHandlerType & dof,
1589 const typename std::vector<typename DoFHandlerType::active_cell_iterator>
1593 DoFHandlerType::space_dimension> *p =
1595 DoFHandlerType::dimension,
1596 DoFHandlerType::space_dimension
> *>(&dof.get_triangulation()))
1602 AssertDimension(cells.size(), dof.get_triangulation().n_active_cells());
1605 const auto n_owned_dofs = dof.n_locally_owned_dofs();
1616 std::vector<bool> already_sorted(n_owned_dofs,
false);
1617 std::vector<types::global_dof_index> cell_dofs;
1619 const auto &owned_dofs = dof.locally_owned_dofs();
1621 unsigned int index = 0;
1623 for (
const auto &cell : cells)
1627 const unsigned int n_cell_dofs = cell->get_fe().n_dofs_per_cell();
1628 cell_dofs.resize(n_cell_dofs);
1630 cell->get_active_or_mg_dof_indices(cell_dofs);
1634 std::sort(cell_dofs.begin(), cell_dofs.end());
1636 for (
const auto dof : cell_dofs)
1638 const auto local_dof = owned_dofs.index_within_set(dof);
1640 !already_sorted[local_dof])
1642 already_sorted[local_dof] =
true;
1643 reverse[index++] = local_dof;
1647 Assert(index == n_owned_dofs,
1649 "Traversing over the given set of cells did not cover all "
1650 "degrees of freedom in the DoFHandler. Does the set of cells "
1651 "not include all active cells?"));
1654 new_indices[reverse[i]] = owned_dofs.nth_index_in_set(i);
1659 template <
typename DoFHandlerType>
1662 DoFHandlerType & dof,
1663 const unsigned int level,
1664 const typename std::vector<typename DoFHandlerType::level_cell_iterator>
1670 std::vector<types::global_dof_index> renumbering(dof.n_dofs(
level));
1671 std::vector<types::global_dof_index> reverse(dof.n_dofs(
level));
1674 dof.renumber_dofs(
level, renumbering);
1679 template <
typename DoFHandlerType>
1682 std::vector<types::global_dof_index> &new_order,
1683 std::vector<types::global_dof_index> &reverse,
1684 const DoFHandlerType & dof,
1685 const unsigned int level,
1686 const typename std::vector<typename DoFHandlerType::level_cell_iterator>
1689 Assert(cells.size() == dof.get_triangulation().n_cells(
level),
1691 dof.get_triangulation().n_cells(
level)));
1697 unsigned int n_global_dofs = dof.n_dofs(
level);
1698 unsigned int n_cell_dofs = dof.get_fe().n_dofs_per_cell();
1700 std::vector<bool> already_sorted(n_global_dofs,
false);
1701 std::vector<types::global_dof_index> cell_dofs(n_cell_dofs);
1705 for (
const auto &cell : cells)
1709 cell->get_active_or_mg_dof_indices(cell_dofs);
1710 std::sort(cell_dofs.begin(), cell_dofs.end());
1712 for (
unsigned int i = 0; i < n_cell_dofs; ++i)
1714 if (!already_sorted[cell_dofs[i]])
1716 already_sorted[cell_dofs[i]] =
true;
1723 "Traversing over the given set of cells did not cover all "
1724 "degrees of freedom in the DoFHandler. Does the set of cells "
1725 "not include all cells of the specified level?"));
1728 new_order[reverse[i]] = i;
1733 template <
typename DoFHandlerType>
1737 const bool dof_wise_renumbering)
1739 std::vector<types::global_dof_index> renumbering(dof.n_dofs());
1740 std::vector<types::global_dof_index> reverse(dof.n_dofs());
1742 renumbering, reverse, dof, direction, dof_wise_renumbering);
1744 dof.renumber_dofs(renumbering);
1749 template <
typename DoFHandlerType>
1752 std::vector<types::global_dof_index> & new_indices,
1753 std::vector<types::global_dof_index> & reverse,
1754 const DoFHandlerType & dof,
1756 const bool dof_wise_renumbering)
1761 DoFHandlerType::space_dimension
> *>(
1762 &dof.get_triangulation()) ==
nullptr),
1765 if (dof_wise_renumbering ==
false)
1767 std::vector<typename DoFHandlerType::active_cell_iterator>
1769 ordered_cells.reserve(dof.get_triangulation().n_active_cells());
1771 DoFHandlerType::space_dimension>
1772 comparator(direction);
1774 typename DoFHandlerType::active_cell_iterator p = dof.begin_active();
1775 typename DoFHandlerType::active_cell_iterator
end = dof.end();
1779 ordered_cells.push_back(p);
1782 std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
1794 const unsigned int n_dofs = dof.n_dofs();
1796 std::pair<Point<DoFHandlerType::space_dimension>,
unsigned int>>
1797 support_point_list(n_dofs);
1800 dof.get_fe_collection();
1801 Assert(fe_collection[0].has_support_points(),
1803 DoFHandlerType::dimension>::ExcFEHasNoSupportPoints());
1805 for (
unsigned int comp = 0; comp < fe_collection.
size(); ++comp)
1807 Assert(fe_collection[comp].has_support_points(),
1809 DoFHandlerType::dimension>::ExcFEHasNoSupportPoints());
1812 fe_collection[comp].get_unit_support_points()));
1815 hp_fe_values(fe_collection,
1816 quadrature_collection,
1819 std::vector<bool> already_touched(n_dofs,
false);
1821 std::vector<types::global_dof_index> local_dof_indices;
1822 typename DoFHandlerType::active_cell_iterator
begin =
1824 typename DoFHandlerType::active_cell_iterator
end = dof.end();
1827 const unsigned int dofs_per_cell =
begin->get_fe().dofs_per_cell;
1828 local_dof_indices.resize(dofs_per_cell);
1832 begin->get_active_or_mg_dof_indices(local_dof_indices);
1833 const std::vector<Point<DoFHandlerType::space_dimension>> &points =
1834 fe_values.get_quadrature_points();
1835 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1836 if (!already_touched[local_dof_indices[i]])
1838 support_point_list[local_dof_indices[i]].first = points[i];
1839 support_point_list[local_dof_indices[i]].second =
1840 local_dof_indices[i];
1841 already_touched[local_dof_indices[i]] =
true;
1847 std::sort(support_point_list.begin(),
1848 support_point_list.end(),
1851 new_indices[support_point_list[i].
second] = i;
1857 template <
typename DoFHandlerType>
1860 const unsigned int level,
1862 const bool dof_wise_renumbering)
1864 std::vector<types::global_dof_index> renumbering(dof.n_dofs(
level));
1865 std::vector<types::global_dof_index> reverse(dof.n_dofs(
level));
1867 renumbering, reverse, dof,
level, direction, dof_wise_renumbering);
1869 dof.renumber_dofs(
level, renumbering);
1874 template <
typename DoFHandlerType>
1877 std::vector<types::global_dof_index> & new_indices,
1878 std::vector<types::global_dof_index> & reverse,
1879 const DoFHandlerType & dof,
1880 const unsigned int level,
1882 const bool dof_wise_renumbering)
1884 if (dof_wise_renumbering ==
false)
1886 std::vector<typename DoFHandlerType::level_cell_iterator> ordered_cells;
1887 ordered_cells.reserve(dof.get_triangulation().n_cells(
level));
1889 DoFHandlerType::space_dimension>
1890 comparator(direction);
1892 typename DoFHandlerType::level_cell_iterator p = dof.begin(
level);
1893 typename DoFHandlerType::level_cell_iterator
end = dof.end(
level);
1897 ordered_cells.push_back(p);
1900 std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
1906 Assert(dof.get_fe().has_support_points(),
1908 DoFHandlerType::dimension>::ExcFEHasNoSupportPoints());
1909 const unsigned int n_dofs = dof.n_dofs(
level);
1911 std::pair<Point<DoFHandlerType::space_dimension>,
unsigned int>>
1912 support_point_list(n_dofs);
1915 dof.get_fe().get_unit_support_points());
1919 std::vector<bool> already_touched(dof.n_dofs(),
false);
1921 const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
1922 std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1923 typename DoFHandlerType::level_cell_iterator
begin = dof.begin(
level);
1924 typename DoFHandlerType::level_cell_iterator
end = dof.end(
level);
1928 DoFHandlerType::dimension,
1929 DoFHandlerType::space_dimension>::cell_iterator &begin_tria =
1931 begin->get_active_or_mg_dof_indices(local_dof_indices);
1932 fe_values.
reinit(begin_tria);
1933 const std::vector<Point<DoFHandlerType::space_dimension>> &points =
1934 fe_values.get_quadrature_points();
1935 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1936 if (!already_touched[local_dof_indices[i]])
1938 support_point_list[local_dof_indices[i]].first = points[i];
1939 support_point_list[local_dof_indices[i]].second =
1940 local_dof_indices[i];
1941 already_touched[local_dof_indices[i]] =
true;
1947 std::sort(support_point_list.begin(),
1948 support_point_list.end(),
1951 new_indices[support_point_list[i].
second] = i;
1985 template <
class DHCellIterator>
1987 operator()(
const DHCellIterator &c1,
const DHCellIterator &c2)
const
1991 return compare(c1, c2, std::integral_constant<int, dim>());
1998 template <
class DHCellIterator,
int xdim>
2001 const DHCellIterator &c2,
2002 std::integral_constant<int, xdim>)
const
2008 return (counter ? (s1 > s2) : (s2 > s1));
2016 template <
class DHCellIterator>
2019 const DHCellIterator &,
2020 std::integral_constant<int, 1>)
const
2023 ExcMessage(
"This operation only makes sense for dim>=2."));
2031 template <
typename DoFHandlerType>
2037 std::vector<types::global_dof_index> renumbering(dof.n_dofs());
2040 dof.renumber_dofs(renumbering);
2045 template <
typename DoFHandlerType>
2048 const DoFHandlerType & dof,
2052 std::vector<typename DoFHandlerType::active_cell_iterator> ordered_cells;
2053 ordered_cells.reserve(dof.get_triangulation().n_active_cells());
2057 typename DoFHandlerType::active_cell_iterator p = dof.begin_active();
2058 typename DoFHandlerType::active_cell_iterator
end = dof.end();
2062 ordered_cells.push_back(p);
2065 std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
2067 std::vector<types::global_dof_index> reverse(new_indices.size());
2073 template <
typename DoFHandlerType>
2076 const unsigned int level,
2080 std::vector<typename DoFHandlerType::level_cell_iterator> ordered_cells;
2081 ordered_cells.reserve(dof.get_triangulation().n_active_cells());
2085 typename DoFHandlerType::level_cell_iterator p = dof.begin(
level);
2086 typename DoFHandlerType::level_cell_iterator
end = dof.end(
level);
2090 ordered_cells.push_back(p);
2093 std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
2100 template <
typename DoFHandlerType>
2104 std::vector<types::global_dof_index> renumbering(
2108 dof_handler.renumber_dofs(renumbering);
2113 template <
typename DoFHandlerType>
2120 std::vector<types::global_dof_index> renumbering(
2121 dof_handler.locally_owned_mg_dofs(
level).n_elements(),
2126 dof_handler.renumber_dofs(
level, renumbering);
2131 template <
typename DoFHandlerType>
2134 const DoFHandlerType & dof_handler)
2137 Assert(new_indices.size() == n_dofs,
2140 std::iota(new_indices.begin(),
2148 ::boost::mt19937 random_number_generator;
2149 for (
unsigned int i = 1; i < n_dofs; ++i)
2152 const unsigned int j =
2153 ::boost::random::uniform_int_distribution<>(0, i)(
2154 random_number_generator);
2158 std::swap(new_indices[i], new_indices[j]);
2164 template <
typename DoFHandlerType>
2167 const DoFHandlerType & dof_handler,
2168 const unsigned int level)
2171 Assert(new_indices.size() == n_dofs,
2174 std::iota(new_indices.begin(),
2182 ::boost::mt19937 random_number_generator;
2183 for (
unsigned int i = 1; i < n_dofs; ++i)
2186 const unsigned int j =
2187 ::boost::random::uniform_int_distribution<>(0, i)(
2188 random_number_generator);
2192 std::swap(new_indices[i], new_indices[j]);
2198 template <
typename DoFHandlerType>
2205 DoFHandlerType::space_dimension
> *>(
2206 &dof_handler.get_triangulation())),
2208 "Parallel triangulations are already enumerated according to their MPI process id."));
2210 std::vector<types::global_dof_index> renumbering(
2214 dof_handler.renumber_dofs(renumbering);
2219 template <
typename DoFHandlerType>
2222 const DoFHandlerType & dof_handler)
2225 Assert(new_dof_indices.size() == n_dofs,
2231 std::vector<types::subdomain_id> subdomain_association(n_dofs);
2233 const unsigned int n_subdomains =
2234 *std::max_element(subdomain_association.begin(),
2235 subdomain_association.end()) +
2245 std::fill(new_dof_indices.begin(),
2246 new_dof_indices.end(),
2252 if (subdomain_association[i] == subdomain)
2256 new_dof_indices[i] = next_free_index;
2262 Assert(std::find(new_dof_indices.begin(),
2263 new_dof_indices.end(),
2273 #include "dof_renumbering.inst"