46#define BOOST_BIND_GLOBAL_PLACEHOLDERS
47#include <boost/config.hpp>
48#include <boost/graph/adjacency_list.hpp>
49#include <boost/graph/bandwidth.hpp>
50#include <boost/graph/cuthill_mckee_ordering.hpp>
51#include <boost/graph/king_ordering.hpp>
52#include <boost/graph/minimum_degree_ordering.hpp>
53#include <boost/graph/properties.hpp>
54#include <boost/random.hpp>
55#include <boost/random/uniform_int_distribution.hpp>
56#undef BOOST_BIND_GLOBAL_PLACEHOLDERS
74 using namespace ::
boost;
77 using Graph = adjacency_list<vecS,
80 property<vertex_color_t,
82 property<vertex_degree_t, int>>>;
83 using Vertex = graph_traits<Graph>::vertex_descriptor;
84 using size_type = graph_traits<Graph>::vertices_size_type;
86 using Pair = std::pair<size_type, size_type>;
92 template <
int dim,
int spacedim>
95 const bool use_constraints,
98 boosttypes::vertex_degree_t>::type
113 for (
unsigned int row = 0; row < dsp.
n_rows(); ++row)
114 for (
unsigned int col = 0; col < dsp.
row_length(row); ++col)
118 boosttypes::graph_traits<boosttypes::Graph>::vertex_iterator ui, ui_end;
120 graph_degree = get(::boost::vertex_degree, graph);
121 for (::boost::tie(ui, ui_end) =
vertices(graph); ui != ui_end; ++ui)
122 graph_degree[*ui] = degree(*ui, graph);
127 template <
int dim,
int spacedim>
130 const bool reversed_numbering,
131 const bool use_constraints)
133 std::vector<types::global_dof_index> renumbering(
147 template <
int dim,
int spacedim>
151 const bool reversed_numbering,
152 const bool use_constraints)
156 boosttypes::vertex_degree_t>::type graph_degree;
161 boosttypes::vertex_index_t>::type index_map =
162 get(::boost::vertex_index, graph);
165 std::vector<boosttypes::Vertex> inv_perm(num_vertices(graph));
167 if (reversed_numbering ==
false)
168 ::boost::cuthill_mckee_ordering(graph,
170 get(::boost::vertex_color, graph),
171 make_degree_map(graph));
173 ::boost::cuthill_mckee_ordering(graph,
175 get(::boost::vertex_color, graph),
176 make_degree_map(graph));
179 new_dof_indices[index_map[inv_perm[c]]] = c;
181 Assert(std::find(new_dof_indices.begin(),
182 new_dof_indices.end(),
189 template <
int dim,
int spacedim>
192 const bool reversed_numbering,
193 const bool use_constraints)
195 std::vector<types::global_dof_index> renumbering(
209 template <
int dim,
int spacedim>
213 const bool reversed_numbering,
214 const bool use_constraints)
218 boosttypes::vertex_degree_t>::type graph_degree;
223 boosttypes::vertex_index_t>::type index_map =
224 get(::boost::vertex_index, graph);
227 std::vector<boosttypes::Vertex> inv_perm(num_vertices(graph));
229 if (reversed_numbering ==
false)
235 new_dof_indices[index_map[inv_perm[c]]] = c;
237 Assert(std::find(new_dof_indices.begin(),
238 new_dof_indices.end(),
245 template <
int dim,
int spacedim>
248 const bool reversed_numbering,
249 const bool use_constraints)
251 std::vector<types::global_dof_index> renumbering(
265 template <
int dim,
int spacedim>
268 std::vector<types::global_dof_index> &new_dof_indices,
270 const bool reversed_numbering,
271 const bool use_constraints)
273 (void)use_constraints;
282 using namespace ::
boost;
287 using Graph = adjacency_list<vecS, vecS, directedS>;
289 int n = dof_handler.
n_dofs();
293 std::vector<::types::global_dof_index> dofs_on_this_cell;
297 const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
299 dofs_on_this_cell.resize(dofs_per_cell);
301 cell->get_active_or_mg_dof_indices(dofs_on_this_cell);
302 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
303 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
304 if (dofs_on_this_cell[i] > dofs_on_this_cell[j])
306 add_edge(dofs_on_this_cell[i], dofs_on_this_cell[j], G);
307 add_edge(dofs_on_this_cell[j], dofs_on_this_cell[i], G);
312 using Vector = std::vector<int>;
315 Vector inverse_perm(n, 0);
320 Vector supernode_sizes(n, 1);
323 ::boost::property_map<Graph, vertex_index_t>::type
id =
324 get(vertex_index, G);
330 minimum_degree_ordering(
332 make_iterator_property_map(degree.
begin(),
id, degree[0]),
335 make_iterator_property_map(supernode_sizes.
begin(),
342 for (
int i = 0; i < n; ++i)
352 if (reversed_numbering ==
true)
357 new_dof_indices.begin());
364 template <
int dim,
int spacedim>
367 const bool reversed_numbering,
368 const bool use_constraints,
369 const std::vector<types::global_dof_index> &starting_indices)
371 std::vector<types::global_dof_index> renumbering(
388 template <
int dim,
int spacedim>
391 std::vector<types::global_dof_index> & new_indices,
393 const bool reversed_numbering,
394 const bool use_constraints,
395 const std::vector<types::global_dof_index> &starting_indices,
396 const unsigned int level)
398 const bool reorder_level_dofs =
415 if (reorder_level_dofs ==
false)
418 locally_relevant_dofs);
427 locally_relevant_dofs);
438 constraints.
reinit(locally_relevant_dofs);
449 locally_owned_dofs.
size());
450 if (reorder_level_dofs ==
false)
462 if (reversed_numbering)
476 if (reorder_level_dofs ==
false)
479 locally_active_dofs);
488 bool needs_locally_active =
false;
489 for (
const auto starting_index : starting_indices)
491 if ((needs_locally_active ==
493 (locally_owned_dofs.
is_element(starting_index) ==
false))
496 locally_active_dofs.
is_element(starting_index),
498 "You specified global degree of freedom " +
500 " as a starting index, but this index is not among the "
501 "locally active ones on this processor, as required "
502 "for this function."));
503 needs_locally_active =
true;
508 (needs_locally_active ? locally_active_dofs : locally_owned_dofs);
519 index_set_to_use.
size(),
521 if (reorder_level_dofs ==
false)
532 std::vector<types::global_dof_index> row_entries;
533 for (
unsigned int i = 0; i < index_set_to_use.
n_elements(); ++i)
537 const unsigned int row_length = dsp.
row_length(row);
539 for (
unsigned int j = 0; j < row_length; ++j)
542 if (col != row && index_set_to_use.
is_element(col))
552 std::vector<types::global_dof_index> local_starting_indices(
553 starting_indices.size());
554 for (
unsigned int i = 0; i < starting_indices.size(); ++i)
555 local_starting_indices[i] =
560 std::vector<types::global_dof_index> my_new_indices(
564 local_starting_indices);
565 if (reversed_numbering)
574 if (needs_locally_active ==
true)
577 IndexSet active_but_not_owned_dofs = locally_active_dofs;
578 active_but_not_owned_dofs.
subtract_set(locally_owned_dofs);
580 std::set<types::global_dof_index> erase_these_indices;
581 for (
const auto p : active_but_not_owned_dofs)
586 erase_these_indices.insert(my_new_indices[index]);
589 Assert(erase_these_indices.size() ==
592 Assert(
static_cast<unsigned int>(
593 std::count(my_new_indices.begin(),
594 my_new_indices.end(),
600 std::vector<types::global_dof_index> translate_indices(
601 my_new_indices.size());
603 std::set<types::global_dof_index>::const_iterator
604 next_erased_index = erase_these_indices.begin();
606 for (
unsigned int i = 0; i < translate_indices.size(); ++i)
607 if ((next_erased_index != erase_these_indices.end()) &&
608 (*next_erased_index == i))
615 translate_indices[i] = next_new_index;
625 new_indices.reserve(locally_owned_dofs.
n_elements());
626 for (
const auto &p : my_new_indices)
631 new_indices.push_back(translate_indices[p]);
637 new_indices = std::move(my_new_indices);
650 template <
int dim,
int spacedim>
653 const unsigned int level,
654 const bool reversed_numbering,
655 const std::vector<types::global_dof_index> &starting_indices)
660 std::vector<types::global_dof_index> new_indices(
679 template <
int dim,
int spacedim>
682 const std::vector<unsigned int> &component_order_arg)
684 std::vector<types::global_dof_index> renumbering(
688 compute_component_wise<dim, spacedim>(renumbering,
705 (result <= dof_handler.
n_dofs())),
713 template <
int dim,
int spacedim>
716 const unsigned int level,
717 const std::vector<unsigned int> &component_order_arg)
722 std::vector<types::global_dof_index> renumbering(
732 compute_component_wise<dim, spacedim>(
733 renumbering, start,
end, component_order_arg,
true);
740 if (renumbering.size() != 0)
746 template <
int dim,
int spacedim,
typename CellIterator>
749 const CellIterator & start,
751 const std::vector<unsigned int> &component_order_arg,
752 const bool is_level_operation)
755 start->get_dof_handler().get_fe_collection());
761 new_indices.resize(0);
768 const IndexSet &locally_owned_dofs =
770 start->get_dof_handler().locally_owned_mg_dofs(start->level()) :
771 start->get_dof_handler().locally_owned_dofs();
775 std::vector<unsigned int> component_order(component_order_arg);
780 if (component_order.size() == 0)
781 for (
unsigned int i = 0; i < fe_collection.
n_components(); ++i)
782 component_order.push_back(i);
788 for (
const unsigned int component : component_order)
796 std::vector<types::global_dof_index> local_dof_indices;
808 std::vector<std::vector<unsigned int>> component_list(fe_collection.
size());
809 for (
unsigned int f = 0; f < fe_collection.
size(); ++f)
813 component_list[f].resize(dofs_per_cell);
814 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
816 component_list[f][i] =
820 const unsigned int comp =
826 component_list[f][i] = component_order[comp];
841 std::vector<std::vector<types::global_dof_index>> component_to_dof_map(
843 for (CellIterator cell = start; cell !=
end; ++cell)
845 if (is_level_operation)
848 if (!cell->is_locally_owned_on_level())
855 if (!cell->is_locally_owned())
859 if (is_level_operation)
861 cell->level() == start->level(),
863 "Multigrid renumbering in compute_component_wise() needs to be applied to a single level!"));
868 const unsigned int fe_index = cell->active_fe_index();
869 const unsigned int dofs_per_cell =
870 fe_collection[fe_index].n_dofs_per_cell();
871 local_dof_indices.resize(dofs_per_cell);
872 cell->get_active_or_mg_dof_indices(local_dof_indices);
874 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
875 if (locally_owned_dofs.
is_element(local_dof_indices[i]))
876 component_to_dof_map[component_list[fe_index][i]].
push_back(
877 local_dof_indices[i]);
903 for (
unsigned int component = 0; component < fe_collection.
n_components();
906 std::sort(component_to_dof_map[component].
begin(),
907 component_to_dof_map[component].
end());
908 component_to_dof_map[component].erase(
909 std::unique(component_to_dof_map[component].
begin(),
910 component_to_dof_map[component].
end()),
911 component_to_dof_map[component].
end());
916 const unsigned int n_buckets = fe_collection.
n_components();
917 std::vector<types::global_dof_index> shifts(n_buckets);
921 &start->get_dof_handler().get_triangulation())))
923#ifdef DEAL_II_WITH_MPI
924 std::vector<types::global_dof_index> local_dof_count(n_buckets);
926 for (
unsigned int c = 0; c < n_buckets; ++c)
927 local_dof_count[c] = component_to_dof_map[c].size();
929 std::vector<types::global_dof_index> prefix_dof_count(n_buckets);
930 const int ierr = MPI_Exscan(local_dof_count.data(),
931 prefix_dof_count.data(),
935 tria->get_communicator());
938 std::vector<types::global_dof_index> global_dof_count(n_buckets);
940 tria->get_communicator(),
945 for (
unsigned int c = 0; c < n_buckets; ++c)
947 shifts[c] = prefix_dof_count[c] + cumulated;
948 cumulated += global_dof_count[c];
958 for (
unsigned int c = 1; c < fe_collection.
n_components(); ++c)
959 shifts[c] = shifts[c - 1] + component_to_dof_map[c - 1].size();
968 for (
unsigned int component = 0; component < fe_collection.
n_components();
971 const typename std::vector<types::global_dof_index>::const_iterator
972 begin_of_component = component_to_dof_map[component].begin(),
973 end_of_component = component_to_dof_map[component].end();
975 next_free_index = shifts[component];
977 for (
typename std::vector<types::global_dof_index>::const_iterator
978 dof_index = begin_of_component;
979 dof_index != end_of_component;
990 return next_free_index;
995 template <
int dim,
int spacedim>
999 std::vector<types::global_dof_index> renumbering(
1020 (result <= dof_handler.
n_dofs())),
1028 template <
int dim,
int spacedim>
1035 std::vector<types::global_dof_index> renumbering(
1058 if (renumbering.size() != 0)
1064 template <
int dim,
int spacedim,
class ITERATOR,
class ENDITERATOR>
1067 const ITERATOR & start,
1068 const ENDITERATOR &
end,
1069 const bool is_level_operation)
1072 start->get_dof_handler().get_fe_collection());
1078 new_indices.resize(0);
1085 const IndexSet &locally_owned_dofs =
1086 is_level_operation ?
1087 start->get_dof_handler().locally_owned_mg_dofs(start->level()) :
1088 start->get_dof_handler().locally_owned_dofs();
1092 std::vector<types::global_dof_index> local_dof_indices;
1097 std::vector<std::vector<types::global_dof_index>> block_list(
1098 fe_collection.
size());
1099 for (
unsigned int f = 0; f < fe_collection.
size(); ++f)
1118 std::vector<std::vector<types::global_dof_index>> block_to_dof_map(
1120 for (ITERATOR cell = start; cell !=
end; ++cell)
1122 if (is_level_operation)
1125 if (!cell->is_locally_owned_on_level())
1132 if (!cell->is_locally_owned())
1136 if (is_level_operation)
1138 cell->level() == start->level(),
1140 "Multigrid renumbering in compute_block_wise() needs to be applied to a single level!"));
1145 const unsigned int fe_index = cell->active_fe_index();
1146 const unsigned int dofs_per_cell =
1147 fe_collection[fe_index].n_dofs_per_cell();
1148 local_dof_indices.resize(dofs_per_cell);
1149 cell->get_active_or_mg_dof_indices(local_dof_indices);
1151 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1152 if (locally_owned_dofs.
is_element(local_dof_indices[i]))
1153 block_to_dof_map[block_list[fe_index][i]].
push_back(
1154 local_dof_indices[i]);
1169 for (
unsigned int block = 0; block < fe_collection.
n_blocks(); ++block)
1171 std::sort(block_to_dof_map[block].
begin(),
1172 block_to_dof_map[block].
end());
1173 block_to_dof_map[block].erase(
1174 std::unique(block_to_dof_map[block].
begin(),
1175 block_to_dof_map[block].
end()),
1176 block_to_dof_map[block].
end());
1181 const unsigned int n_buckets = fe_collection.
n_blocks();
1182 std::vector<types::global_dof_index> shifts(n_buckets);
1186 &start->get_dof_handler().get_triangulation())))
1188#ifdef DEAL_II_WITH_MPI
1189 std::vector<types::global_dof_index> local_dof_count(n_buckets);
1191 for (
unsigned int c = 0; c < n_buckets; ++c)
1192 local_dof_count[c] = block_to_dof_map[c].size();
1194 std::vector<types::global_dof_index> prefix_dof_count(n_buckets);
1195 const int ierr = MPI_Exscan(local_dof_count.data(),
1196 prefix_dof_count.data(),
1200 tria->get_communicator());
1203 std::vector<types::global_dof_index> global_dof_count(n_buckets);
1205 tria->get_communicator(),
1210 for (
unsigned int c = 0; c < n_buckets; ++c)
1212 shifts[c] = prefix_dof_count[c] + cumulated;
1213 cumulated += global_dof_count[c];
1223 for (
unsigned int c = 1; c < fe_collection.
n_blocks(); ++c)
1224 shifts[c] = shifts[c - 1] + block_to_dof_map[c - 1].size();
1233 for (
unsigned int block = 0; block < fe_collection.
n_blocks(); ++block)
1235 const typename std::vector<types::global_dof_index>::const_iterator
1236 begin_of_component = block_to_dof_map[block].begin(),
1237 end_of_component = block_to_dof_map[block].end();
1239 next_free_index = shifts[block];
1241 for (
typename std::vector<types::global_dof_index>::const_iterator
1242 dof_index = begin_of_component;
1243 dof_index != end_of_component;
1254 return next_free_index;
1266 template <
int dim,
class CellIteratorType>
1268 compute_hierarchical_recursive(
1271 const CellIteratorType & cell,
1272 const IndexSet & locally_owned_dof_indices,
1273 std::vector<types::global_dof_index> &new_indices)
1276 next_free_dof_offset;
1278 if (cell->has_children())
1281 for (
unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell;
1283 current_next_free_dof_offset =
1284 compute_hierarchical_recursive<dim>(current_next_free_dof_offset,
1287 locally_owned_dof_indices,
1307 if (cell->is_locally_owned())
1310 const unsigned int dofs_per_cell =
1311 cell->get_fe().n_dofs_per_cell();
1312 std::vector<types::global_dof_index> local_dof_indices(
1314 cell->get_dof_indices(local_dof_indices);
1333 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1334 if (locally_owned_dof_indices.
is_element(local_dof_indices[i]))
1338 const unsigned int idx =
1340 local_dof_indices[i]);
1344 my_starting_index + current_next_free_dof_offset;
1345 ++current_next_free_dof_offset;
1351 return current_next_free_dof_offset;
1357 template <
int dim,
int spacedim>
1361 std::vector<types::global_dof_index> renumbering(
1384#ifdef DEAL_II_WITH_MPI
1387 const int ierr = MPI_Exscan(&locally_owned_size,
1392 tria->get_communicator());
1401#ifdef DEAL_II_WITH_P4EST
1406 for (
unsigned int c = 0; c < tria->n_cells(0); ++c)
1408 const unsigned int coarse_cell_index =
1409 tria->get_p4est_tree_to_coarse_cell_permutation()[c];
1412 this_cell(tria, 0, coarse_cell_index, &dof_handler);
1414 next_free_dof_offset =
1415 compute_hierarchical_recursive<dim>(next_free_dof_offset,
1430 dof_handler.
begin(0);
1431 cell != dof_handler.
end(0);
1433 next_free_dof_offset =
1434 compute_hierarchical_recursive<dim>(next_free_dof_offset,
1447 (next_free_dof_offset <= dof_handler.
n_dofs())),
1451 Assert(std::find(renumbering.begin(),
1461 template <
int dim,
int spacedim>
1464 const std::vector<bool> & selected_dofs)
1466 std::vector<types::global_dof_index> renumbering(
1475 template <
int dim,
int spacedim>
1478 const std::vector<bool> & selected_dofs,
1479 const unsigned int level)
1484 std::vector<types::global_dof_index> renumbering(
1496 template <
int dim,
int spacedim>
1499 std::vector<types::global_dof_index> &new_indices,
1501 const std::vector<bool> & selected_dofs)
1504 Assert(selected_dofs.size() == n_dofs,
1509 Assert(new_indices.size() == n_dofs,
1513 std::count(selected_dofs.begin(), selected_dofs.end(),
false);
1518 if (selected_dofs[i] ==
false)
1520 new_indices[i] = next_unselected;
1525 new_indices[i] = next_selected;
1534 template <
int dim,
int spacedim>
1537 std::vector<types::global_dof_index> &new_indices,
1539 const std::vector<bool> & selected_dofs,
1540 const unsigned int level)
1545 const unsigned int n_dofs = dof_handler.
n_dofs(
level);
1546 Assert(selected_dofs.size() == n_dofs,
1551 Assert(new_indices.size() == n_dofs,
1554 const unsigned int n_selected_dofs =
1555 std::count(selected_dofs.begin(), selected_dofs.end(),
false);
1557 unsigned int next_unselected = 0;
1558 unsigned int next_selected = n_selected_dofs;
1559 for (
unsigned int i = 0; i < n_dofs; ++i)
1560 if (selected_dofs[i] ==
false)
1562 new_indices[i] = next_unselected;
1567 new_indices[i] = next_selected;
1576 template <
int dim,
int spacedim>
1583 std::vector<types::global_dof_index> renumbering(
1592 template <
int dim,
int spacedim>
1595 std::vector<types::global_dof_index> &new_indices,
1596 std::vector<types::global_dof_index> &reverse,
1598 const typename std::vector<
1623 std::vector<bool> already_sorted(n_owned_dofs,
false);
1624 std::vector<types::global_dof_index> cell_dofs;
1628 unsigned int index = 0;
1630 for (
const auto &cell : cells)
1634 const unsigned int n_cell_dofs = cell->get_fe().n_dofs_per_cell();
1635 cell_dofs.resize(n_cell_dofs);
1637 cell->get_active_or_mg_dof_indices(cell_dofs);
1641 std::sort(cell_dofs.begin(), cell_dofs.end());
1643 for (
const auto dof : cell_dofs)
1645 const auto local_dof = owned_dofs.index_within_set(dof);
1647 !already_sorted[local_dof])
1649 already_sorted[local_dof] =
true;
1650 reverse[index++] = local_dof;
1654 Assert(index == n_owned_dofs,
1656 "Traversing over the given set of cells did not cover all "
1657 "degrees of freedom in the DoFHandler. Does the set of cells "
1658 "not include all active cells?"));
1661 new_indices[reverse[i]] = owned_dofs.nth_index_in_set(i);
1666 template <
int dim,
int spacedim>
1669 const unsigned int level,
1670 const typename std::vector<
1676 std::vector<types::global_dof_index> renumbering(dof.
n_dofs(
level));
1677 std::vector<types::global_dof_index> reverse(dof.
n_dofs(
level));
1685 template <
int dim,
int spacedim>
1688 std::vector<types::global_dof_index> &new_order,
1689 std::vector<types::global_dof_index> &reverse,
1691 const unsigned int level,
1692 const typename std::vector<
1706 std::vector<bool> already_sorted(n_global_dofs,
false);
1707 std::vector<types::global_dof_index> cell_dofs(n_cell_dofs);
1711 for (
const auto &cell : cells)
1715 cell->get_active_or_mg_dof_indices(cell_dofs);
1716 std::sort(cell_dofs.begin(), cell_dofs.end());
1718 for (
unsigned int i = 0; i < n_cell_dofs; ++i)
1720 if (!already_sorted[cell_dofs[i]])
1722 already_sorted[cell_dofs[i]] =
true;
1729 "Traversing over the given set of cells did not cover all "
1730 "degrees of freedom in the DoFHandler. Does the set of cells "
1731 "not include all cells of the specified level?"));
1734 new_order[reverse[i]] = i;
1739 template <
int dim,
int spacedim>
1743 const bool dof_wise_renumbering)
1745 std::vector<types::global_dof_index> renumbering(dof.
n_dofs());
1746 std::vector<types::global_dof_index> reverse(dof.
n_dofs());
1748 renumbering, reverse, dof, direction, dof_wise_renumbering);
1755 template <
int dim,
int spacedim>
1758 std::vector<types::global_dof_index> &reverse,
1761 const bool dof_wise_renumbering)
1767 if (dof_wise_renumbering ==
false)
1769 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
1775 comparator(direction);
1778 ordered_cells.push_back(cell);
1780 std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
1792 const unsigned int n_dofs = dof.
n_dofs();
1793 std::vector<std::pair<Point<spacedim>,
unsigned int>>
1794 support_point_list(n_dofs);
1797 Assert(fe_collection[0].has_support_points(),
1800 for (
unsigned int comp = 0; comp < fe_collection.
size(); ++comp)
1802 Assert(fe_collection[comp].has_support_points(),
1808 quadrature_collection,
1811 std::vector<bool> already_touched(n_dofs,
false);
1813 std::vector<types::global_dof_index> local_dof_indices;
1817 const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
1818 local_dof_indices.resize(dofs_per_cell);
1819 hp_fe_values.
reinit(cell);
1822 cell->get_active_or_mg_dof_indices(local_dof_indices);
1823 const std::vector<Point<spacedim>> &points =
1824 fe_values.get_quadrature_points();
1825 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1826 if (!already_touched[local_dof_indices[i]])
1828 support_point_list[local_dof_indices[i]].first = points[i];
1829 support_point_list[local_dof_indices[i]].second =
1830 local_dof_indices[i];
1831 already_touched[local_dof_indices[i]] =
true;
1836 std::sort(support_point_list.begin(),
1837 support_point_list.end(),
1840 new_indices[support_point_list[i].
second] = i;
1846 template <
int dim,
int spacedim>
1849 const unsigned int level,
1851 const bool dof_wise_renumbering)
1853 std::vector<types::global_dof_index> renumbering(dof.
n_dofs(
level));
1854 std::vector<types::global_dof_index> reverse(dof.
n_dofs(
level));
1856 renumbering, reverse, dof,
level, direction, dof_wise_renumbering);
1863 template <
int dim,
int spacedim>
1866 std::vector<types::global_dof_index> &reverse,
1868 const unsigned int level,
1870 const bool dof_wise_renumbering)
1872 if (dof_wise_renumbering ==
false)
1874 std::vector<typename DoFHandler<dim, spacedim>::level_cell_iterator>
1880 comparator(direction);
1889 ordered_cells.push_back(p);
1892 std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
1901 std::vector<std::pair<Point<spacedim>,
unsigned int>>
1902 support_point_list(n_dofs);
1909 std::vector<bool> already_touched(dof.
n_dofs(),
false);
1912 std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1920 &begin_tria =
begin;
1921 begin->get_active_or_mg_dof_indices(local_dof_indices);
1922 fe_values.
reinit(begin_tria);
1923 const std::vector<Point<spacedim>> &points =
1925 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1926 if (!already_touched[local_dof_indices[i]])
1928 support_point_list[local_dof_indices[i]].first = points[i];
1929 support_point_list[local_dof_indices[i]].second =
1930 local_dof_indices[i];
1931 already_touched[local_dof_indices[i]] =
true;
1936 std::sort(support_point_list.begin(),
1937 support_point_list.end(),
1940 new_indices[support_point_list[i].
second] = i;
1974 template <
class DHCellIterator>
1976 operator()(
const DHCellIterator &c1,
const DHCellIterator &c2)
const
1980 return compare(c1, c2, std::integral_constant<int, dim>());
1987 template <
class DHCellIterator,
int xdim>
1990 const DHCellIterator &c2,
1991 std::integral_constant<int, xdim>)
const
1997 return (
counter ? (s1 > s2) : (s2 > s1));
2005 template <
class DHCellIterator>
2008 const DHCellIterator &,
2009 std::integral_constant<int, 1>)
const
2012 ExcMessage(
"This operation only makes sense for dim>=2."));
2020 template <
int dim,
int spacedim>
2026 std::vector<types::global_dof_index> renumbering(dof.
n_dofs());
2034 template <
int dim,
int spacedim>
2041 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
2047 ordered_cells.push_back(cell);
2049 std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
2051 std::vector<types::global_dof_index> reverse(new_indices.size());
2057 template <
int dim,
int spacedim>
2060 const unsigned int level,
2064 std::vector<typename DoFHandler<dim, spacedim>::level_cell_iterator>
2076 ordered_cells.push_back(p);
2079 std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
2086 template <
int dim,
int spacedim>
2090 std::vector<types::global_dof_index> renumbering(
2099 template <
int dim,
int spacedim>
2106 std::vector<types::global_dof_index> renumbering(
2117 template <
int dim,
int spacedim>
2123 Assert(new_indices.size() == n_dofs,
2126 std::iota(new_indices.begin(),
2134 ::boost::mt19937 random_number_generator;
2135 for (
unsigned int i = 1; i < n_dofs; ++i)
2138 const unsigned int j =
2139 ::boost::random::uniform_int_distribution<>(0, i)(
2140 random_number_generator);
2144 std::swap(new_indices[i], new_indices[j]);
2150 template <
int dim,
int spacedim>
2154 const unsigned int level)
2157 Assert(new_indices.size() == n_dofs,
2160 std::iota(new_indices.begin(),
2168 ::boost::mt19937 random_number_generator;
2169 for (
unsigned int i = 1; i < n_dofs; ++i)
2172 const unsigned int j =
2173 ::boost::random::uniform_int_distribution<>(0, i)(
2174 random_number_generator);
2178 std::swap(new_indices[i], new_indices[j]);
2184 template <
int dim,
int spacedim>
2192 "Parallel triangulations are already enumerated according to their MPI process id."));
2194 std::vector<types::global_dof_index> renumbering(
2203 template <
int dim,
int spacedim>
2209 Assert(new_dof_indices.size() == n_dofs,
2215 std::vector<types::subdomain_id> subdomain_association(n_dofs);
2217 const unsigned int n_subdomains =
2218 *std::max_element(subdomain_association.begin(),
2219 subdomain_association.end()) +
2229 std::fill(new_dof_indices.begin(),
2230 new_dof_indices.end(),
2236 if (subdomain_association[i] == subdomain)
2240 new_dof_indices[i] = next_free_index;
2246 Assert(std::find(new_dof_indices.begin(),
2247 new_dof_indices.end(),
2257#include "dof_renumbering.inst"
void reinit(const IndexSet &local_constraints=IndexSet())
unsigned int first_selected_component(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
cell_iterator end() const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
void renumber_dofs(const std::vector< types::global_dof_index > &new_numbers)
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
const Triangulation< dim, spacedim > & get_triangulation() const
const IndexSet & locally_owned_dofs() const
active_cell_iterator begin_active(const unsigned int level=0) const
types::global_dof_index n_dofs() const
cell_iterator begin(const unsigned int level=0) const
types::global_dof_index n_locally_owned_dofs() const
const std::vector< Point< spacedim > > & get_quadrature_points() const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
unsigned int n_dofs_per_cell() const
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
const ComponentMask & get_nonzero_components(const unsigned int i) const
bool has_support_points() const
const std::vector< Point< dim > > & get_unit_support_points() const
bool is_primitive() const
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
size_type index_within_set(const size_type global_index) const
size_type n_elements() const
bool is_element(const size_type index) const
void subtract_set(const IndexSet &other)
size_type nth_index_in_set(const size_type local_index) const
unsigned int n_active_cells() const
unsigned int n_cells() const
unsigned int size() const
void push_back(const FiniteElement< dim, spacedim > &new_fe)
unsigned int n_blocks() const
unsigned int n_components() const
const FEValuesType & get_present_fe_values() const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, lda > > &cell, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
void push_back(const Quadrature< dim_in > &new_quadrature)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
@ update_quadrature_points
Transformed quadrature points.
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
std::string to_string(const T &t)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDoFHandlerNotInitialized()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
void add_entries(const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_unique_and_sorted=false)
size_type row_length(const size_type row) const
size_type column_number(const size_type row, const size_type index) const
void make_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity_pattern, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
Expression atan2(const Expression &y, const Expression &x)
graph_traits< Graph >::vertices_size_type size_type
graph_traits< Graph >::vertex_descriptor Vertex
adjacency_list< vecS, vecS, undirectedS, property< vertex_color_t, default_color_type, property< vertex_degree_t, int > > > Graph
std::pair< size_type, size_type > Pair
void create_graph(const DoFHandler< dim, spacedim > &dof_handler, const bool use_constraints, boosttypes::Graph &graph, boosttypes::property_map< boosttypes::Graph, boosttypes::vertex_degree_t >::type &graph_degree)
void compute_Cuthill_McKee(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandler< dim, spacedim > &, const bool reversed_numbering=false, const bool use_constraints=false)
void Cuthill_McKee(DoFHandler< dim, spacedim > &dof_handler, const bool reversed_numbering=false, const bool use_constraints=false)
void king_ordering(DoFHandler< dim, spacedim > &dof_handler, const bool reversed_numbering=false, const bool use_constraints=false)
void compute_king_ordering(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandler< dim, spacedim > &, const bool reversed_numbering=false, const bool use_constraints=false)
void compute_minimum_degree(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandler< dim, spacedim > &, const bool reversed_numbering=false, const bool use_constraints=false)
void minimum_degree(DoFHandler< dim, spacedim > &dof_handler, const bool reversed_numbering=false, const bool use_constraints=false)
void compute_subdomain_wise(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandler< dim, spacedim > &dof_handler)
void subdomain_wise(DoFHandler< dim, spacedim > &dof_handler)
void hierarchical(DoFHandler< dim, spacedim > &dof_handler)
void compute_random(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandler< dim, spacedim > &dof_handler)
void component_wise(DoFHandler< dim, spacedim > &dof_handler, const std::vector< unsigned int > &target_component=std::vector< unsigned int >())
void downstream(DoFHandler< dim, spacedim > &dof_handler, const Tensor< 1, spacedim > &direction, const bool dof_wise_renumbering=false)
types::global_dof_index compute_component_wise(std::vector< types::global_dof_index > &new_dof_indices, const CellIterator &start, const typename identity< CellIterator >::type &end, const std::vector< unsigned int > &target_component, const bool is_level_operation)
void compute_Cuthill_McKee(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandler< dim, spacedim > &, const bool reversed_numbering=false, const bool use_constraints=false, const std::vector< types::global_dof_index > &starting_indices=std::vector< types::global_dof_index >(), const unsigned int level=numbers::invalid_unsigned_int)
void block_wise(DoFHandler< dim, spacedim > &dof_handler)
void Cuthill_McKee(DoFHandler< dim, spacedim > &dof_handler, const bool reversed_numbering=false, const bool use_constraints=false, const std::vector< types::global_dof_index > &starting_indices=std::vector< types::global_dof_index >())
void compute_sort_selected_dofs_back(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandler< dim, spacedim > &dof_handler, const std::vector< bool > &selected_dofs)
void sort_selected_dofs_back(DoFHandler< dim, spacedim > &dof_handler, const std::vector< bool > &selected_dofs)
void compute_cell_wise(std::vector< types::global_dof_index > &renumbering, std::vector< types::global_dof_index > &inverse_renumbering, const DoFHandler< dim, spacedim > &dof_handler, const std::vector< typename DoFHandler< dim, spacedim >::active_cell_iterator > &cell_order)
void clockwise_dg(DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim > ¢er, const bool counter=false)
void random(DoFHandler< dim, spacedim > &dof_handler)
void compute_downstream(std::vector< types::global_dof_index > &new_dof_indices, std::vector< types::global_dof_index > &reverse, const DoFHandler< dim, spacedim > &dof_handler, const Tensor< 1, spacedim > &direction, const bool dof_wise_renumbering)
void cell_wise(DoFHandler< dim, spacedim > &dof_handler, const std::vector< typename DoFHandler< dim, spacedim >::active_cell_iterator > &cell_order)
types::global_dof_index compute_block_wise(std::vector< types::global_dof_index > &new_dof_indices, const ITERATOR &start, const ENDITERATOR &end, bool is_level_operation)
void compute_clockwise_dg(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim > ¢er, const bool counter)
void swap(MemorySpaceData< Number, MemorySpace > &, MemorySpaceData< Number, MemorySpace > &)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
T sum(const T &t, const MPI_Comm &mpi_communicator)
std::vector< Integer > reverse_permutation(const std::vector< Integer > &permutation)
void copy(const T *begin, const T *end, U *dest)
static const unsigned int invalid_unsigned_int
const types::global_dof_index invalid_dof_index
bool compare(const DHCellIterator &c1, const DHCellIterator &c2, std::integral_constant< int, xdim >) const
bool compare(const DHCellIterator &, const DHCellIterator &, std::integral_constant< int, 1 >) const
ClockCells(const Point< dim > ¢er, bool counter)
const Point< dim > & center
bool operator()(const DHCellIterator &c1, const DHCellIterator &c2) const
#define DEAL_II_DOF_INDEX_MPI_TYPE