49#define BOOST_BIND_GLOBAL_PLACEHOLDERS
50#include <boost/config.hpp>
51#include <boost/graph/adjacency_list.hpp>
52#include <boost/graph/bandwidth.hpp>
53#include <boost/graph/cuthill_mckee_ordering.hpp>
54#include <boost/graph/king_ordering.hpp>
55#include <boost/graph/minimum_degree_ordering.hpp>
56#include <boost/graph/properties.hpp>
57#include <boost/random.hpp>
58#include <boost/random/uniform_int_distribution.hpp>
60#undef BOOST_BIND_GLOBAL_PLACEHOLDERS
78 using namespace ::
boost;
80 using Graph = adjacency_list<vecS,
83 property<vertex_color_t,
85 property<vertex_degree_t, int>>>;
86 using Vertex = graph_traits<Graph>::vertex_descriptor;
87 using size_type = graph_traits<Graph>::vertices_size_type;
89 using Pair = std::pair<size_type, size_type>;
95 template <
int dim,
int spacedim>
98 const bool use_constraints,
101 boosttypes::vertex_degree_t>::type
116 for (
unsigned int row = 0; row < dsp.
n_rows(); ++row)
117 for (
unsigned int col = 0; col < dsp.
row_length(row); ++col)
121 boosttypes::graph_traits<boosttypes::Graph>::vertex_iterator ui, ui_end;
123 graph_degree = get(::boost::vertex_degree, graph);
124 for (::boost::tie(ui, ui_end) =
vertices(graph); ui != ui_end; ++ui)
125 graph_degree[*ui] = degree(*ui, graph);
130 template <
int dim,
int spacedim>
133 const bool reversed_numbering,
134 const bool use_constraints)
136 std::vector<types::global_dof_index> renumbering(
150 template <
int dim,
int spacedim>
154 const bool reversed_numbering,
155 const bool use_constraints)
159 boosttypes::vertex_degree_t>::type graph_degree;
164 boosttypes::vertex_index_t>::type index_map =
165 get(::boost::vertex_index, graph);
168 std::vector<boosttypes::Vertex> inv_perm(num_vertices(graph));
170 if (reversed_numbering ==
false)
171 ::boost::cuthill_mckee_ordering(graph,
173 get(::boost::vertex_color, graph),
174 make_degree_map(graph));
176 ::boost::cuthill_mckee_ordering(graph,
178 get(::boost::vertex_color, graph),
179 make_degree_map(graph));
182 new_dof_indices[index_map[inv_perm[c]]] = c;
184 Assert(std::find(new_dof_indices.begin(),
185 new_dof_indices.end(),
192 template <
int dim,
int spacedim>
195 const bool reversed_numbering,
196 const bool use_constraints)
198 std::vector<types::global_dof_index> renumbering(
212 template <
int dim,
int spacedim>
216 const bool reversed_numbering,
217 const bool use_constraints)
221 boosttypes::vertex_degree_t>::type graph_degree;
226 boosttypes::vertex_index_t>::type index_map =
227 get(::boost::vertex_index, graph);
230 std::vector<boosttypes::Vertex> inv_perm(num_vertices(graph));
232 if (reversed_numbering ==
false)
235 ::boost::king_ordering(graph, inv_perm.begin());
238 new_dof_indices[index_map[inv_perm[c]]] = c;
240 Assert(std::find(new_dof_indices.begin(),
241 new_dof_indices.end(),
248 template <
int dim,
int spacedim>
251 const bool reversed_numbering,
252 const bool use_constraints)
254 std::vector<types::global_dof_index> renumbering(
268 template <
int dim,
int spacedim>
271 std::vector<types::global_dof_index> &new_dof_indices,
273 const bool reversed_numbering,
274 const bool use_constraints)
276 (void)use_constraints;
285 using namespace ::
boost;
290 using Graph = adjacency_list<vecS, vecS, directedS>;
292 int n = dof_handler.
n_dofs();
296 std::vector<::types::global_dof_index> dofs_on_this_cell;
300 const unsigned int dofs_per_cell = cell->get_fe().n_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)
355 if (reversed_numbering ==
true)
356 std::copy(perm.
begin(), perm.
end(), new_dof_indices.begin());
358 std::copy(inverse_perm.
begin(),
360 new_dof_indices.begin());
367 template <
int dim,
int spacedim>
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(
391 template <
int dim,
int spacedim>
394 std::vector<types::global_dof_index> & new_indices,
396 const bool reversed_numbering,
397 const bool use_constraints,
398 const std::vector<types::global_dof_index> &starting_indices,
399 const unsigned int level)
401 const bool reorder_level_dofs =
418 if (reorder_level_dofs ==
false)
421 locally_relevant_dofs);
430 locally_relevant_dofs);
441 constraints.
reinit(locally_relevant_dofs);
452 locally_owned_dofs.
size());
453 if (reorder_level_dofs ==
false)
465 if (reversed_numbering)
479 if (reorder_level_dofs ==
false)
482 locally_active_dofs);
491 bool needs_locally_active =
false;
492 for (
const auto starting_index : starting_indices)
494 if ((needs_locally_active ==
496 (locally_owned_dofs.
is_element(starting_index) ==
false))
499 locally_active_dofs.
is_element(starting_index),
501 "You specified global degree of freedom " +
502 std::to_string(starting_index) +
503 " as a starting index, but this index is not among the "
504 "locally active ones on this processor, as required "
505 "for this function."));
506 needs_locally_active =
true;
511 (needs_locally_active ? locally_active_dofs : locally_owned_dofs);
522 index_set_to_use.
size(),
524 if (reorder_level_dofs ==
false)
535 std::vector<types::global_dof_index> row_entries;
536 for (
unsigned int i = 0; i < index_set_to_use.
n_elements(); ++i)
540 const unsigned int row_length = dsp.
row_length(row);
542 for (
unsigned int j = 0; j < row_length; ++j)
545 if (col != row && index_set_to_use.
is_element(col))
555 std::vector<types::global_dof_index> local_starting_indices(
556 starting_indices.size());
557 for (
unsigned int i = 0; i < starting_indices.size(); ++i)
558 local_starting_indices[i] =
563 std::vector<types::global_dof_index> my_new_indices(
567 local_starting_indices);
568 if (reversed_numbering)
577 if (needs_locally_active ==
true)
580 IndexSet active_but_not_owned_dofs = locally_active_dofs;
581 active_but_not_owned_dofs.
subtract_set(locally_owned_dofs);
583 std::set<types::global_dof_index> erase_these_indices;
584 for (
const auto p : active_but_not_owned_dofs)
589 erase_these_indices.insert(my_new_indices[index]);
592 Assert(erase_these_indices.size() ==
595 Assert(
static_cast<unsigned int>(
596 std::count(my_new_indices.begin(),
597 my_new_indices.end(),
603 std::vector<types::global_dof_index> translate_indices(
604 my_new_indices.size());
606 std::set<types::global_dof_index>::const_iterator
607 next_erased_index = erase_these_indices.begin();
609 for (
unsigned int i = 0; i < translate_indices.size(); ++i)
610 if ((next_erased_index != erase_these_indices.end()) &&
611 (*next_erased_index == i))
618 translate_indices[i] = next_new_index;
628 new_indices.reserve(locally_owned_dofs.
n_elements());
629 for (
const auto &p : my_new_indices)
634 new_indices.push_back(translate_indices[p]);
640 new_indices = std::move(my_new_indices);
653 template <
int dim,
int spacedim>
656 const unsigned int level,
657 const bool reversed_numbering,
658 const std::vector<types::global_dof_index> &starting_indices)
663 std::vector<types::global_dof_index> new_indices(
682 template <
int dim,
int spacedim>
685 const std::vector<unsigned int> &component_order_arg)
687 std::vector<types::global_dof_index> renumbering(
691 compute_component_wise<dim, spacedim>(renumbering,
713 (result <= dof_handler.
n_dofs())),
721 template <
int dim,
int spacedim>
724 const unsigned int level,
725 const std::vector<unsigned int> &component_order_arg)
730 std::vector<types::global_dof_index> renumbering(
740 compute_component_wise<dim, spacedim>(
741 renumbering, start, end, component_order_arg,
true);
754 template <
int dim,
int spacedim,
typename CellIterator>
757 const CellIterator & start,
759 const std::vector<unsigned int> &component_order_arg,
760 const bool is_level_operation)
763 start->get_dof_handler().get_fe_collection());
769 new_indices.resize(0);
776 const IndexSet &locally_owned_dofs =
778 start->get_dof_handler().locally_owned_mg_dofs(start->level()) :
779 start->get_dof_handler().locally_owned_dofs();
783 std::vector<unsigned int> component_order(component_order_arg);
788 if (component_order.size() == 0)
789 for (
unsigned int i = 0; i < fe_collection.
n_components(); ++i)
790 component_order.push_back(i);
796 for (
const unsigned int component : component_order)
804 std::vector<types::global_dof_index> local_dof_indices;
816 std::vector<std::vector<unsigned int>> component_list(fe_collection.
size());
817 for (
unsigned int f = 0; f < fe_collection.
size(); ++f)
821 component_list[f].resize(dofs_per_cell);
822 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
824 component_list[f][i] =
828 const unsigned int comp =
834 component_list[f][i] = component_order[comp];
849 std::vector<std::vector<types::global_dof_index>> component_to_dof_map(
851 for (CellIterator cell = start; cell != end; ++cell)
853 if (is_level_operation)
856 if (!cell->is_locally_owned_on_level())
863 if (!cell->is_locally_owned())
867 if (is_level_operation)
869 cell->level() == start->level(),
871 "Multigrid renumbering in compute_component_wise() needs to be applied to a single level!"));
876 const unsigned int fe_index = cell->active_fe_index();
877 const unsigned int dofs_per_cell =
878 fe_collection[fe_index].n_dofs_per_cell();
879 local_dof_indices.resize(dofs_per_cell);
880 cell->get_active_or_mg_dof_indices(local_dof_indices);
882 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
883 if (locally_owned_dofs.
is_element(local_dof_indices[i]))
884 component_to_dof_map[component_list[fe_index][i]].
push_back(
885 local_dof_indices[i]);
911 for (
unsigned int component = 0; component < fe_collection.
n_components();
914 std::sort(component_to_dof_map[component].begin(),
915 component_to_dof_map[component].end());
916 component_to_dof_map[component].erase(
917 std::unique(component_to_dof_map[component].begin(),
918 component_to_dof_map[component].end()),
919 component_to_dof_map[component].end());
924 const unsigned int n_buckets = fe_collection.
n_components();
925 std::vector<types::global_dof_index> shifts(n_buckets);
929 &start->get_dof_handler().get_triangulation())))
931#ifdef DEAL_II_WITH_MPI
932 std::vector<types::global_dof_index> local_dof_count(n_buckets);
934 for (
unsigned int c = 0; c < n_buckets; ++c)
935 local_dof_count[c] = component_to_dof_map[c].size();
937 std::vector<types::global_dof_index> prefix_dof_count(n_buckets);
938 const int ierr = MPI_Exscan(local_dof_count.data(),
939 prefix_dof_count.data(),
946 std::vector<types::global_dof_index> global_dof_count(n_buckets);
953 for (
unsigned int c = 0; c < n_buckets; ++c)
955 shifts[c] = prefix_dof_count[c] + cumulated;
956 cumulated += global_dof_count[c];
966 for (
unsigned int c = 1; c < fe_collection.
n_components(); ++c)
967 shifts[c] = shifts[c - 1] + component_to_dof_map[c - 1].size();
976 for (
unsigned int component = 0; component < fe_collection.
n_components();
979 next_free_index = shifts[component];
982 component_to_dof_map[component])
994 return next_free_index;
999 template <
int dim,
int spacedim>
1003 std::vector<types::global_dof_index> renumbering(
1024 (result <= dof_handler.
n_dofs())),
1032 template <
int dim,
int spacedim>
1039 std::vector<types::global_dof_index> renumbering(
1068 template <
int dim,
int spacedim,
class ITERATOR,
class ENDITERATOR>
1071 const ITERATOR & start,
1072 const ENDITERATOR & end,
1073 const bool is_level_operation)
1076 start->get_dof_handler().get_fe_collection());
1082 new_indices.resize(0);
1089 const IndexSet &locally_owned_dofs =
1090 is_level_operation ?
1091 start->get_dof_handler().locally_owned_mg_dofs(start->level()) :
1092 start->get_dof_handler().locally_owned_dofs();
1096 std::vector<types::global_dof_index> local_dof_indices;
1101 std::vector<std::vector<types::global_dof_index>> block_list(
1102 fe_collection.
size());
1103 for (
unsigned int f = 0; f < fe_collection.
size(); ++f)
1122 std::vector<std::vector<types::global_dof_index>> block_to_dof_map(
1124 for (ITERATOR cell = start; cell != end; ++cell)
1126 if (is_level_operation)
1129 if (!cell->is_locally_owned_on_level())
1136 if (!cell->is_locally_owned())
1140 if (is_level_operation)
1142 cell->level() == start->level(),
1144 "Multigrid renumbering in compute_block_wise() needs to be applied to a single level!"));
1149 const unsigned int fe_index = cell->active_fe_index();
1150 const unsigned int dofs_per_cell =
1151 fe_collection[fe_index].n_dofs_per_cell();
1152 local_dof_indices.resize(dofs_per_cell);
1153 cell->get_active_or_mg_dof_indices(local_dof_indices);
1155 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1156 if (locally_owned_dofs.
is_element(local_dof_indices[i]))
1157 block_to_dof_map[block_list[fe_index][i]].
push_back(
1158 local_dof_indices[i]);
1173 for (
unsigned int block = 0; block < fe_collection.
n_blocks(); ++block)
1175 std::sort(block_to_dof_map[block].begin(),
1176 block_to_dof_map[block].end());
1177 block_to_dof_map[block].erase(
1178 std::unique(block_to_dof_map[block].begin(),
1179 block_to_dof_map[block].end()),
1180 block_to_dof_map[block].end());
1185 const unsigned int n_buckets = fe_collection.
n_blocks();
1186 std::vector<types::global_dof_index> shifts(n_buckets);
1190 &start->get_dof_handler().get_triangulation())))
1192#ifdef DEAL_II_WITH_MPI
1193 std::vector<types::global_dof_index> local_dof_count(n_buckets);
1195 for (
unsigned int c = 0; c < n_buckets; ++c)
1196 local_dof_count[c] = block_to_dof_map[c].size();
1198 std::vector<types::global_dof_index> prefix_dof_count(n_buckets);
1199 const int ierr = MPI_Exscan(local_dof_count.data(),
1200 prefix_dof_count.data(),
1207 std::vector<types::global_dof_index> global_dof_count(n_buckets);
1214 for (
unsigned int c = 0; c < n_buckets; ++c)
1216 shifts[c] = prefix_dof_count[c] + cumulated;
1217 cumulated += global_dof_count[c];
1227 for (
unsigned int c = 1; c < fe_collection.
n_blocks(); ++c)
1228 shifts[c] = shifts[c - 1] + block_to_dof_map[c - 1].size();
1237 for (
unsigned int block = 0; block < fe_collection.
n_blocks(); ++block)
1239 const typename std::vector<types::global_dof_index>::const_iterator
1240 begin_of_component = block_to_dof_map[block].begin(),
1241 end_of_component = block_to_dof_map[block].end();
1243 next_free_index = shifts[block];
1245 for (
typename std::vector<types::global_dof_index>::const_iterator
1246 dof_index = begin_of_component;
1247 dof_index != end_of_component;
1258 return next_free_index;
1270 template <
int dim,
class CellIteratorType>
1272 compute_hierarchical_recursive(
1275 const CellIteratorType & cell,
1276 const IndexSet & locally_owned_dof_indices,
1277 std::vector<types::global_dof_index> &new_indices)
1280 next_free_dof_offset;
1282 if (cell->has_children())
1285 for (
unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell;
1287 current_next_free_dof_offset =
1288 compute_hierarchical_recursive<dim>(current_next_free_dof_offset,
1291 locally_owned_dof_indices,
1311 if (cell->is_locally_owned())
1314 const unsigned int dofs_per_cell =
1315 cell->get_fe().n_dofs_per_cell();
1316 std::vector<types::global_dof_index> local_dof_indices(
1318 cell->get_dof_indices(local_dof_indices);
1337 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1338 if (locally_owned_dof_indices.
is_element(local_dof_indices[i]))
1342 const unsigned int idx =
1344 local_dof_indices[i]);
1348 my_starting_index + current_next_free_dof_offset;
1349 ++current_next_free_dof_offset;
1355 return current_next_free_dof_offset;
1361 template <
int dim,
int spacedim>
1365 std::vector<types::global_dof_index> renumbering(
1388#ifdef DEAL_II_WITH_MPI
1391 const int ierr = MPI_Exscan(&locally_owned_size,
1405#ifdef DEAL_II_WITH_P4EST
1410 for (
unsigned int c = 0; c <
tria->
n_cells(0); ++c)
1412 const unsigned int coarse_cell_index =
1413 tria->get_p4est_tree_to_coarse_cell_permutation()[c];
1416 this_cell(
tria, 0, coarse_cell_index, &dof_handler);
1418 next_free_dof_offset =
1419 compute_hierarchical_recursive<dim>(next_free_dof_offset,
1434 dof_handler.
begin(0);
1435 cell != dof_handler.
end(0);
1437 next_free_dof_offset =
1438 compute_hierarchical_recursive<dim>(next_free_dof_offset,
1451 (next_free_dof_offset <= dof_handler.
n_dofs())),
1455 Assert(std::find(renumbering.begin(),
1465 template <
int dim,
int spacedim>
1468 const std::vector<bool> & selected_dofs)
1470 std::vector<types::global_dof_index> renumbering(
1479 template <
int dim,
int spacedim>
1482 const std::vector<bool> & selected_dofs,
1483 const unsigned int level)
1488 std::vector<types::global_dof_index> renumbering(
1500 template <
int dim,
int spacedim>
1503 std::vector<types::global_dof_index> &new_indices,
1505 const std::vector<bool> & selected_dofs)
1508 Assert(selected_dofs.size() == n_dofs,
1513 Assert(new_indices.size() == n_dofs,
1517 std::count(selected_dofs.begin(), selected_dofs.end(),
false);
1522 if (selected_dofs[i] ==
false)
1524 new_indices[i] = next_unselected;
1529 new_indices[i] = next_selected;
1538 template <
int dim,
int spacedim>
1541 std::vector<types::global_dof_index> &new_indices,
1543 const std::vector<bool> & selected_dofs,
1544 const unsigned int level)
1549 const unsigned int n_dofs = dof_handler.
n_dofs(
level);
1550 Assert(selected_dofs.size() == n_dofs,
1555 Assert(new_indices.size() == n_dofs,
1558 const unsigned int n_selected_dofs =
1559 std::count(selected_dofs.begin(), selected_dofs.end(),
false);
1561 unsigned int next_unselected = 0;
1562 unsigned int next_selected = n_selected_dofs;
1563 for (
unsigned int i = 0; i < n_dofs; ++i)
1564 if (selected_dofs[i] ==
false)
1566 new_indices[i] = next_unselected;
1571 new_indices[i] = next_selected;
1580 template <
int dim,
int spacedim>
1587 std::vector<types::global_dof_index> renumbering(
1596 template <
int dim,
int spacedim>
1599 std::vector<types::global_dof_index> &new_indices,
1600 std::vector<types::global_dof_index> &reverse,
1602 const typename std::vector<
1627 std::vector<bool> already_sorted(n_owned_dofs,
false);
1628 std::vector<types::global_dof_index> cell_dofs;
1632 unsigned int index = 0;
1634 for (
const auto &cell : cells)
1638 const unsigned int n_cell_dofs = cell->get_fe().n_dofs_per_cell();
1639 cell_dofs.resize(n_cell_dofs);
1641 cell->get_active_or_mg_dof_indices(cell_dofs);
1645 std::sort(cell_dofs.begin(), cell_dofs.end());
1647 for (
const auto dof : cell_dofs)
1649 const auto local_dof = owned_dofs.index_within_set(dof);
1651 !already_sorted[local_dof])
1653 already_sorted[local_dof] =
true;
1654 reverse[index++] = local_dof;
1658 Assert(index == n_owned_dofs,
1660 "Traversing over the given set of cells did not cover all "
1661 "degrees of freedom in the DoFHandler. Does the set of cells "
1662 "not include all active cells?"));
1665 new_indices[reverse[i]] = owned_dofs.nth_index_in_set(i);
1670 template <
int dim,
int spacedim>
1673 const unsigned int level,
1674 const typename std::vector<
1680 std::vector<types::global_dof_index> renumbering(dof.
n_dofs(
level));
1681 std::vector<types::global_dof_index> reverse(dof.
n_dofs(
level));
1689 template <
int dim,
int spacedim>
1692 std::vector<types::global_dof_index> &new_order,
1693 std::vector<types::global_dof_index> &reverse,
1695 const unsigned int level,
1696 const typename std::vector<
1710 std::vector<bool> already_sorted(n_global_dofs,
false);
1711 std::vector<types::global_dof_index> cell_dofs(n_cell_dofs);
1713 unsigned int global_index = 0;
1715 for (
const auto &cell : cells)
1719 cell->get_active_or_mg_dof_indices(cell_dofs);
1720 std::sort(cell_dofs.begin(), cell_dofs.end());
1722 for (
unsigned int i = 0; i < n_cell_dofs; ++i)
1724 if (!already_sorted[cell_dofs[i]])
1726 already_sorted[cell_dofs[i]] =
true;
1727 reverse[global_index++] = cell_dofs[i];
1731 Assert(global_index == n_global_dofs,
1733 "Traversing over the given set of cells did not cover all "
1734 "degrees of freedom in the DoFHandler. Does the set of cells "
1735 "not include all cells of the specified level?"));
1738 new_order[reverse[i]] = i;
1743 template <
int dim,
int spacedim>
1747 const bool dof_wise_renumbering)
1749 std::vector<types::global_dof_index> renumbering(dof.
n_dofs());
1750 std::vector<types::global_dof_index> reverse(dof.
n_dofs());
1752 renumbering, reverse, dof, direction, dof_wise_renumbering);
1759 template <
int dim,
int spacedim>
1762 std::vector<types::global_dof_index> &reverse,
1765 const bool dof_wise_renumbering)
1771 if (dof_wise_renumbering ==
false)
1773 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
1779 comparator(direction);
1782 ordered_cells.push_back(cell);
1784 std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
1796 const unsigned int n_dofs = dof.
n_dofs();
1797 std::vector<std::pair<Point<spacedim>,
unsigned int>>
1798 support_point_list(n_dofs);
1801 Assert(fe_collection[0].has_support_points(),
1804 for (
unsigned int comp = 0; comp < fe_collection.
size(); ++comp)
1806 Assert(fe_collection[comp].has_support_points(),
1812 quadrature_collection,
1815 std::vector<bool> already_touched(n_dofs,
false);
1817 std::vector<types::global_dof_index> local_dof_indices;
1821 const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
1822 local_dof_indices.resize(dofs_per_cell);
1823 hp_fe_values.
reinit(cell);
1826 cell->get_active_or_mg_dof_indices(local_dof_indices);
1827 const std::vector<Point<spacedim>> &points =
1829 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1830 if (!already_touched[local_dof_indices[i]])
1832 support_point_list[local_dof_indices[i]].first = points[i];
1833 support_point_list[local_dof_indices[i]].second =
1834 local_dof_indices[i];
1835 already_touched[local_dof_indices[i]] =
true;
1840 std::sort(support_point_list.begin(),
1841 support_point_list.end(),
1844 new_indices[support_point_list[i].
second] = i;
1850 template <
int dim,
int spacedim>
1853 const unsigned int level,
1855 const bool dof_wise_renumbering)
1857 std::vector<types::global_dof_index> renumbering(dof.
n_dofs(
level));
1858 std::vector<types::global_dof_index> reverse(dof.
n_dofs(
level));
1860 renumbering, reverse, dof,
level, direction, dof_wise_renumbering);
1867 template <
int dim,
int spacedim>
1870 std::vector<types::global_dof_index> &reverse,
1872 const unsigned int level,
1874 const bool dof_wise_renumbering)
1876 if (dof_wise_renumbering ==
false)
1878 std::vector<typename DoFHandler<dim, spacedim>::level_cell_iterator>
1884 comparator(direction);
1893 ordered_cells.push_back(p);
1896 std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
1905 std::vector<std::pair<Point<spacedim>,
unsigned int>>
1906 support_point_list(n_dofs);
1913 std::vector<bool> already_touched(dof.
n_dofs(),
false);
1916 std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1921 for (; begin != end; ++begin)
1924 &begin_tria = begin;
1925 begin->get_active_or_mg_dof_indices(local_dof_indices);
1926 fe_values.
reinit(begin_tria);
1927 const std::vector<Point<spacedim>> &points =
1929 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1930 if (!already_touched[local_dof_indices[i]])
1932 support_point_list[local_dof_indices[i]].first = points[i];
1933 support_point_list[local_dof_indices[i]].second =
1934 local_dof_indices[i];
1935 already_touched[local_dof_indices[i]] =
true;
1940 std::sort(support_point_list.begin(),
1941 support_point_list.end(),
1944 new_indices[support_point_list[i].
second] = i;
1978 template <
class DHCellIterator>
1980 operator()(
const DHCellIterator &c1,
const DHCellIterator &c2)
const
1984 return compare(c1, c2, std::integral_constant<int, dim>());
1991 template <
class DHCellIterator,
int xdim>
1994 const DHCellIterator &c2,
1995 std::integral_constant<int, xdim>)
const
1999 const double s1 = std::atan2(
v1[0],
v1[1]);
2000 const double s2 = std::atan2(v2[0], v2[1]);
2001 return (
counter ? (s1 > s2) : (s2 > s1));
2009 template <
class DHCellIterator>
2012 const DHCellIterator &,
2013 std::integral_constant<int, 1>)
const
2016 ExcMessage(
"This operation only makes sense for dim>=2."));
2024 template <
int dim,
int spacedim>
2030 std::vector<types::global_dof_index> renumbering(dof.
n_dofs());
2038 template <
int dim,
int spacedim>
2045 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
2051 ordered_cells.push_back(cell);
2053 std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
2055 std::vector<types::global_dof_index> reverse(new_indices.size());
2061 template <
int dim,
int spacedim>
2064 const unsigned int level,
2068 std::vector<typename DoFHandler<dim, spacedim>::level_cell_iterator>
2080 ordered_cells.push_back(p);
2083 std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
2090 template <
int dim,
int spacedim>
2094 std::vector<types::global_dof_index> renumbering(
2103 template <
int dim,
int spacedim>
2110 std::vector<types::global_dof_index> renumbering(
2121 template <
int dim,
int spacedim>
2127 Assert(new_indices.size() == n_dofs,
2130 std::iota(new_indices.begin(),
2138 ::boost::mt19937 random_number_generator;
2139 for (
unsigned int i = 1; i < n_dofs; ++i)
2142 const unsigned int j =
2143 ::boost::random::uniform_int_distribution<>(0, i)(
2144 random_number_generator);
2148 std::swap(new_indices[i], new_indices[j]);
2154 template <
int dim,
int spacedim>
2158 const unsigned int level)
2161 Assert(new_indices.size() == n_dofs,
2164 std::iota(new_indices.begin(),
2172 ::boost::mt19937 random_number_generator;
2173 for (
unsigned int i = 1; i < n_dofs; ++i)
2176 const unsigned int j =
2177 ::boost::random::uniform_int_distribution<>(0, i)(
2178 random_number_generator);
2182 std::swap(new_indices[i], new_indices[j]);
2188 template <
int dim,
int spacedim>
2196 "Parallel triangulations are already enumerated according to their MPI process id."));
2198 std::vector<types::global_dof_index> renumbering(
2207 template <
int dim,
int spacedim>
2213 Assert(new_dof_indices.size() == n_dofs,
2219 std::vector<types::subdomain_id> subdomain_association(n_dofs);
2221 const unsigned int n_subdomains =
2222 *std::max_element(subdomain_association.begin(),
2223 subdomain_association.end()) +
2233 std::fill(new_dof_indices.begin(),
2234 new_dof_indices.end(),
2240 if (subdomain_association[i] == subdomain)
2244 new_dof_indices[i] = next_free_index;
2250 Assert(std::find(new_dof_indices.begin(),
2251 new_dof_indices.end(),
2258 template <
int dim,
int spacedim>
2262 std::vector<types::global_dof_index> renumbering(
2275 template <
int dim,
int spacedim>
2278 std::vector<types::global_dof_index> &new_dof_indices,
2282 Assert(new_dof_indices.size() == n_dofs,
2294 std::vector<types::global_dof_index> component_renumbering(
2296 compute_component_wise<dim, spacedim>(component_renumbering,
2299 std::vector<unsigned int>(),
2304 const std::vector<types::global_dof_index> dofs_per_component =
2306 for (
const auto &dpc : dofs_per_component)
2310 const unsigned int n_components =
2314 dof_handler.
n_dofs() / n_components;
2320 if (component_renumbering.size() == 0)
2322 new_dof_indices.resize(0);
2325 std::fill(new_dof_indices.begin(),
2326 new_dof_indices.end(),
2334 std::vector<types::global_dof_index> component_dofs(
2335 local_dofs_per_component);
2336 for (
unsigned int component = 0; component < n_components; ++component)
2338 for (std::size_t i = 0; i < local_dofs_per_component; ++i)
2340 component_renumbering[n_components * i + component];
2341 component_renumbered_dofs.
add_indices(component_dofs.begin(),
2342 component_dofs.end());
2344 component_renumbered_dofs.
compress();
2348 component_renumbered_dofs2.
add_indices(component_renumbering.begin(),
2349 component_renumbering.end());
2350 Assert(component_renumbered_dofs2 == component_renumbered_dofs,
2357 AssertThrow(fe.dofs_per_cell == 0 || fe.has_support_points(),
2359 for (
unsigned int i = 0; i < fe.n_base_elements(); ++i)
2361 fe.base_element(0).get_unit_support_points() ==
2362 fe.base_element(i).get_unit_support_points(),
2364 "All base elements should have the same support points."));
2367 std::vector<types::global_dof_index> component_to_nodal(
2371 std::vector<types::global_dof_index> cell_dofs;
2372 std::vector<types::global_dof_index> component_renumbered_cell_dofs;
2376 auto next_dof_it = locally_owned_dofs.
begin();
2378 if (cell->is_locally_owned())
2383 cell->get_dof_indices(cell_dofs);
2390 if (locally_owned_dofs.
is_element(cell_dofs[i]))
2392 const auto local_index =
2394 component_renumbered_cell_dofs[i] =
2395 component_renumbering[local_index];
2399 component_renumbered_cell_dofs[i] =
2408 component_renumbered_cell_dofs[i]))
2410 for (
unsigned int component = 0;
2417 const auto local_index =
2419 component_renumbered_cell_dofs[i] +
2420 dofs_per_component * component);
2422 if (component_to_nodal[local_index] ==
2425 component_to_nodal[local_index] = *next_dof_it;
2436 const auto local_index =
2438 new_dof_indices[i] = component_to_nodal[local_index];
2447 typename VectorizedArrayType>
2453 const std::vector<types::global_dof_index> new_global_numbers =
2463 template <
int dim,
int spacedim,
typename Number,
typename AdditionalDataType>
2467 const AdditionalDataType & matrix_free_data)
2469 const std::vector<types::global_dof_index> new_global_numbers =
2476 dof_handler.
renumber_dofs(matrix_free_data.mg_level, new_global_numbers);
2481 template <
int dim,
int spacedim,
typename Number,
typename AdditionalDataType>
2482 std::vector<types::global_dof_index>
2486 const AdditionalDataType & matrix_free_data)
2488 AdditionalDataType my_mf_data = matrix_free_data;
2489 my_mf_data.initialize_mapping =
false;
2490 my_mf_data.tasks_parallel_scheme = AdditionalDataType::none;
2492 typename AdditionalDataType::MatrixFreeType separate_matrix_free;
2511 std::vector<std::vector<unsigned int>>
2512 group_dofs_by_rank_access(
2513 const ::Utilities::MPI::Partitioner &partitioner)
2517 std::vector<unsigned int> touch_count(partitioner.locally_owned_size());
2518 for (
const auto &p : partitioner.import_indices())
2519 for (
unsigned int i = p.first; i < p.second; ++i)
2523 std::vector<std::vector<unsigned int>> result(1);
2524 for (
unsigned int i = 0; i < touch_count.size(); ++i)
2525 if (touch_count[i] == 0)
2526 result.back().push_back(i);
2531 std::map<unsigned int, std::vector<unsigned int>>
2532 multiple_ranks_access_dof;
2533 const std::vector<std::pair<unsigned int, unsigned int>> &import_targets =
2534 partitioner.import_targets();
2535 auto it = partitioner.import_indices().begin();
2536 for (
const std::pair<unsigned int, unsigned int> &proc : import_targets)
2538 result.emplace_back();
2539 unsigned int count_dofs = 0;
2540 while (count_dofs < proc.second)
2542 for (
unsigned int i = it->first; i < it->
second;
2545 if (touch_count[i] == 1)
2546 result.back().push_back(i);
2548 multiple_ranks_access_dof[i].push_back(proc.first);
2558 std::map<std::vector<unsigned int>,
2559 std::vector<unsigned int>,
2560 std::function<
bool(
const std::vector<unsigned int> &,
2561 const std::vector<unsigned int> &)>>
2562 dofs_by_rank{[](
const std::vector<unsigned int> &a,
2563 const std::vector<unsigned int> &
b) {
2564 if (a.size() <
b.size())
2566 if (a.size() ==
b.size())
2568 for (
unsigned int i = 0; i < a.size(); ++i)
2571 else if (a[i] > b[i])
2576 for (
const auto &entry : multiple_ranks_access_dof)
2577 dofs_by_rank[entry.second].push_back(entry.first);
2579 for (
const auto &procs : dofs_by_rank)
2580 result.push_back(procs.second);
2591 template <
int dim,
typename Number,
typename VectorizedArrayType>
2592 std::pair<std::vector<unsigned int>, std::vector<unsigned char>>
2593 compute_mf_numbering(
2595 const unsigned int component)
2599 const unsigned int n_comp =
2602 matrix_free.
get_dof_handler(component).get_fe().n_base_elements() == 1,
2610 const unsigned int fe_degree =
2612 const unsigned int nn = fe_degree - 1;
2621 std::array<std::pair<unsigned int, unsigned int>,
2626 dofs_on_objects[0] = std::make_pair(0U, 1U);
2627 dofs_on_objects[1] = std::make_pair(2 * n_comp, nn);
2628 dofs_on_objects[2] = std::make_pair(n_comp, 1U);
2632 dofs_on_objects[0] = std::make_pair(0U, 1U);
2633 dofs_on_objects[1] = std::make_pair(n_comp * (4 + 2 * nn), nn);
2634 dofs_on_objects[2] = std::make_pair(n_comp, 1U);
2635 dofs_on_objects[3] = std::make_pair(n_comp * 4, nn);
2636 dofs_on_objects[4] = std::make_pair(n_comp * (4 + 4 * nn), nn * nn);
2637 dofs_on_objects[5] = std::make_pair(n_comp * (4 + 1 * nn), nn);
2638 dofs_on_objects[6] = std::make_pair(2 * n_comp, 1U);
2639 dofs_on_objects[7] = std::make_pair(n_comp * (4 + 3 * nn), nn);
2640 dofs_on_objects[8] = std::make_pair(3 * n_comp, 1U);
2644 dofs_on_objects[0] = std::make_pair(0U, 1U);
2645 dofs_on_objects[1] = std::make_pair(n_comp * (8 + 2 * nn), nn);
2646 dofs_on_objects[2] = std::make_pair(n_comp, 1U);
2647 dofs_on_objects[3] = std::make_pair(n_comp * 8, nn);
2648 dofs_on_objects[4] =
2649 std::make_pair(n_comp * (8 + 12 * nn + 4 * nn * nn), nn * nn);
2650 dofs_on_objects[5] = std::make_pair(n_comp * (8 + 1 * nn), nn);
2651 dofs_on_objects[6] = std::make_pair(n_comp * 2, 1U);
2652 dofs_on_objects[7] = std::make_pair(n_comp * (8 + 3 * nn), nn);
2653 dofs_on_objects[8] = std::make_pair(n_comp * 3, 1U);
2654 dofs_on_objects[9] = std::make_pair(n_comp * (8 + 8 * nn), nn);
2655 dofs_on_objects[10] =
2656 std::make_pair(n_comp * (8 + 12 * nn + 2 * nn * nn), nn * nn);
2657 dofs_on_objects[11] = std::make_pair(n_comp * (8 + 9 * nn), nn);
2658 dofs_on_objects[12] = std::make_pair(n_comp * (8 + 12 * nn), nn * nn);
2659 dofs_on_objects[13] =
2660 std::make_pair(n_comp * (8 + 12 * nn + 6 * nn * nn), nn * nn * nn);
2661 dofs_on_objects[14] =
2662 std::make_pair(n_comp * (8 + 12 * nn + 1 * nn * nn), nn * nn);
2663 dofs_on_objects[15] = std::make_pair(n_comp * (8 + 10 * nn), nn);
2664 dofs_on_objects[16] =
2665 std::make_pair(n_comp * (8 + 12 * nn + 3 * nn * nn), nn * nn);
2666 dofs_on_objects[17] = std::make_pair(n_comp * (8 + 11 * nn), nn);
2667 dofs_on_objects[18] = std::make_pair(n_comp * 4, 1U);
2668 dofs_on_objects[19] = std::make_pair(n_comp * (8 + 6 * nn), nn);
2669 dofs_on_objects[20] = std::make_pair(n_comp * 5, 1U);
2670 dofs_on_objects[21] = std::make_pair(n_comp * (8 + 4 * nn), nn);
2671 dofs_on_objects[22] =
2672 std::make_pair(n_comp * (8 + 12 * nn + 5 * nn * nn), nn * nn);
2673 dofs_on_objects[23] = std::make_pair(n_comp * (8 + 5 * nn), nn);
2674 dofs_on_objects[24] = std::make_pair(n_comp * 6, 1U);
2675 dofs_on_objects[25] = std::make_pair(n_comp * (8 + 7 * nn), nn);
2676 dofs_on_objects[26] = std::make_pair(n_comp * 7, 1U);
2681 std::vector<unsigned int> & result,
2682 unsigned int &counter_dof_numbers) {
2683 if (owned_dofs.is_element(dof_index))
2685 const unsigned int local_dof_index =
2686 owned_dofs.index_within_set(dof_index);
2689 result[local_dof_index] = counter_dof_numbers++;
2693 unsigned int counter_dof_numbers = 0;
2694 std::vector<unsigned int> local_dofs, dofs_extracted;
2695 std::vector<types::global_dof_index> dof_indices(
2712 std::vector<unsigned int> dof_numbers_mf_order(
2714 std::vector<unsigned char> touch_count(owned_dofs.n_elements());
2716 const auto operation_on_cell_range =
2719 const unsigned int &,
2720 const std::pair<unsigned int, unsigned int> &cell_range) {
2722 for (
unsigned int cell = cell_range.first; cell < cell_range.second;
2726 for (
unsigned int v = 0;
2733 ->get_dof_indices(dof_indices);
2736 ->get_mg_dof_indices(dof_indices);
2738 for (
unsigned int a = 0; a < dofs_on_objects.size(); ++a)
2740 const auto &r = dofs_on_objects[a];
2741 if (a == 10 || a == 16)
2744 for (
unsigned int i1 = 0; i1 < nn; ++i1)
2745 for (
unsigned int i0 = 0; i0 < nn; ++i0)
2746 for (
unsigned int c = 0; c < n_comp; ++c)
2747 renumber_func(dof_indices[r.first + r.second * c +
2750 dof_numbers_mf_order,
2751 counter_dof_numbers);
2753 for (
unsigned int i = 0; i < r.second; ++i)
2754 for (
unsigned int c = 0; c < n_comp; ++c)
2756 dof_indices[r.first + r.second * c + i],
2758 dof_numbers_mf_order,
2759 counter_dof_numbers);
2766 dofs_extracted, cell);
2767 local_dofs.insert(local_dofs.end(),
2768 dofs_extracted.begin(),
2769 dofs_extracted.end());
2772 std::sort(local_dofs.begin(), local_dofs.end());
2773 local_dofs.erase(std::unique(local_dofs.begin(), local_dofs.end()),
2776 for (
unsigned int i : local_dofs)
2777 if (i < touch_count.size())
2785 "version of MatrixFree::cell_loop"));
2787 matrix_free.template cell_loop<unsigned int, unsigned int>(
2788 operation_on_cell_range, counter_dof_numbers, counter_dof_numbers);
2792 return std::make_pair(dof_numbers_mf_order, touch_count);
2802 typename VectorizedArrayType>
2803 std::vector<types::global_dof_index>
2809 ExcMessage(
"You need to set up indices in MatrixFree "
2810 "to be able to compute a renumbering!"));
2813 unsigned int component = 0;
2814 for (; component < matrix_free.
n_components(); ++component)
2819 ExcMessage(
"Could not locate the given DoFHandler in MatrixFree"));
2832 const std::vector<std::vector<unsigned int>> dofs_by_rank_access =
2833 group_dofs_by_rank_access(
2836 const std::pair<std::vector<unsigned int>, std::vector<unsigned char>>
2837 local_numbering = compute_mf_numbering(matrix_free, component);
2842 std::vector<unsigned int> new_numbers;
2843 new_numbers.reserve(owned_dofs.
n_elements());
2849 const std::vector<unsigned int> &purely_local_dofs = dofs_by_rank_access[0];
2850 for (
unsigned int i : purely_local_dofs)
2851 if (local_numbering.second[i] == 1)
2852 new_numbers.push_back(i);
2854 const auto comparator = [&](
const unsigned int a,
const unsigned int b) {
2855 return (local_numbering.first[a] < local_numbering.first[b]);
2858 std::sort(new_numbers.begin(), new_numbers.end(), comparator);
2859 unsigned int sorted_size = new_numbers.size();
2863 for (
auto i : purely_local_dofs)
2864 if (local_numbering.second[i] > 1)
2865 new_numbers.push_back(i);
2866 std::sort(new_numbers.begin() + sorted_size, new_numbers.end(), comparator);
2867 sorted_size = new_numbers.size();
2870 for (
unsigned int chunk = 1; chunk < dofs_by_rank_access.size(); ++chunk)
2871 for (
auto i : dofs_by_rank_access[chunk])
2872 new_numbers.push_back(i);
2873 std::sort(new_numbers.begin() + sorted_size, new_numbers.end(), comparator);
2874 sorted_size = new_numbers.size();
2877 for (
auto i : purely_local_dofs)
2878 if (local_numbering.second[i] == 0)
2879 new_numbers.push_back(i);
2880 std::sort(new_numbers.begin() + sorted_size, new_numbers.end(), comparator);
2884 std::vector<::types::global_dof_index> new_global_numbers(
2886 for (
unsigned int i = 0; i < new_numbers.size(); ++i)
2889 return new_global_numbers;
2897#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
typename LevelSelector::cell_iterator level_cell_iterator
cell_iterator begin(const unsigned int level=0) const
MPI_Comm get_communicator() 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
unsigned int n_components() const
const unsigned int dofs_per_cell
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
ElementIterator begin() const
void subtract_set(const IndexSet &other)
size_type nth_index_in_set(const size_type local_index) const
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
const internal::MatrixFreeFunctions::TaskInfo & get_task_info() const
unsigned int n_active_entries_per_cell_batch(const unsigned int cell_batch_index) const
unsigned int get_mg_level() const
const internal::MatrixFreeFunctions::DoFInfo & get_dof_info(const unsigned int dof_handler_index_component=0) const
const DoFHandler< dim > & get_dof_handler(const unsigned int dof_handler_index=0) const
DoFHandler< dim >::cell_iterator get_cell_iterator(const unsigned int cell_batch_index, const unsigned int lane_index, const unsigned int dof_handler_index=0) const
bool indices_initialized() const
unsigned int n_components() const
virtual MPI_Comm get_communicator() 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)
#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)
#define AssertThrow(cond, exc)
typename ActiveSelector::cell_iterator cell_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
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)
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_support_point_wise(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandler< dim, spacedim > &dof_handler)
void compute_subdomain_wise(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandler< dim, spacedim > &dof_handler)
void matrix_free_data_locality(DoFHandler< dim, spacedim > &dof_handler, const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free)
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 support_point_wise(DoFHandler< dim, spacedim > &dof_handler)
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)
std::vector< types::global_dof_index > compute_matrix_free_data_locality(const DoFHandler< dim, spacedim > &dof_handler, const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free)
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)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
T sum(const T &t, const MPI_Comm &mpi_communicator)
T max(const T &t, const MPI_Comm &mpi_communicator)
constexpr T pow(const T base, const int iexp)
std::vector< Integer > reverse_permutation(const std::vector< Integer > &permutation)
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
void get_dof_indices_on_cell_batch(std::vector< unsigned int > &local_indices, const unsigned int cell_batch, const bool with_constraints=true) const
std::shared_ptr< const Utilities::MPI::Partitioner > vector_partitioner
TasksParallelScheme scheme
const ::Triangulation< dim, spacedim > & tria
#define DEAL_II_DOF_INDEX_MPI_TYPE