Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
dof_renumbering.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
18 #include <deal.II/base/types.h>
19 #include <deal.II/base/utilities.h>
20 
22 
26 #include <deal.II/dofs/dof_tools.h>
27 
28 #include <deal.II/fe/fe.h>
29 
30 #include <deal.II/grid/tria.h>
32 
33 #include <deal.II/hp/dof_handler.h>
35 #include <deal.II/hp/fe_values.h>
36 
41 
43 
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
57 
58 #include <algorithm>
59 #include <cmath>
60 #include <functional>
61 #include <map>
62 #include <vector>
63 
64 
66 
67 namespace DoFRenumbering
68 {
69  namespace boost
70  {
71  namespace boosttypes
72  {
73  using namespace ::boost;
74  using namespace std;
75 
76  using Graph = adjacency_list<vecS,
77  vecS,
78  undirectedS,
79  property<vertex_color_t,
80  default_color_type,
81  property<vertex_degree_t, int>>>;
82  using Vertex = graph_traits<Graph>::vertex_descriptor;
83  using size_type = graph_traits<Graph>::vertices_size_type;
84 
85  using Pair = std::pair<size_type, size_type>;
86  } // namespace boosttypes
87 
88 
89  namespace internal
90  {
91  template <typename DoFHandlerType>
92  void
93  create_graph(const DoFHandlerType &dof_handler,
94  const bool use_constraints,
95  boosttypes::Graph & graph,
96  boosttypes::property_map<boosttypes::Graph,
97  boosttypes::vertex_degree_t>::type
98  &graph_degree)
99  {
100  {
101  // create intermediate sparsity pattern (faster than directly
102  // submitting indices)
103  AffineConstraints<double> constraints;
104  if (use_constraints)
105  DoFTools::make_hanging_node_constraints(dof_handler, constraints);
106  constraints.close();
107  DynamicSparsityPattern dsp(dof_handler.n_dofs(),
108  dof_handler.n_dofs());
109  DoFTools::make_sparsity_pattern(dof_handler, dsp, constraints);
110 
111  // submit the entries to the boost graph
112  for (unsigned int row = 0; row < dsp.n_rows(); ++row)
113  for (unsigned int col = 0; col < dsp.row_length(row); ++col)
114  add_edge(row, dsp.column_number(row, col), graph);
115  }
116 
117  boosttypes::graph_traits<boosttypes::Graph>::vertex_iterator ui, ui_end;
118 
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);
122  }
123  } // namespace internal
124 
125 
126  template <typename DoFHandlerType>
127  void
128  Cuthill_McKee(DoFHandlerType &dof_handler,
129  const bool reversed_numbering,
130  const bool use_constraints)
131  {
132  std::vector<types::global_dof_index> renumbering(
133  dof_handler.n_dofs(), numbers::invalid_dof_index);
134  compute_Cuthill_McKee(renumbering,
135  dof_handler,
136  reversed_numbering,
137  use_constraints);
138 
139  // actually perform renumbering;
140  // this is dimension specific and
141  // thus needs an own function
142  dof_handler.renumber_dofs(renumbering);
143  }
144 
145 
146  template <typename DoFHandlerType>
147  void
148  compute_Cuthill_McKee(std::vector<types::global_dof_index> &new_dof_indices,
149  const DoFHandlerType & dof_handler,
150  const bool reversed_numbering,
151  const bool use_constraints)
152  {
153  boosttypes::Graph graph(dof_handler.n_dofs());
154  boosttypes::property_map<boosttypes::Graph,
155  boosttypes::vertex_degree_t>::type graph_degree;
156 
157  internal::create_graph(dof_handler, use_constraints, graph, graph_degree);
158 
159  boosttypes::property_map<boosttypes::Graph,
160  boosttypes::vertex_index_t>::type index_map =
161  get(::boost::vertex_index, graph);
162 
163 
164  std::vector<boosttypes::Vertex> inv_perm(num_vertices(graph));
165 
166  if (reversed_numbering == false)
167  ::boost::cuthill_mckee_ordering(graph,
168  inv_perm.rbegin(),
169  get(::boost::vertex_color, graph),
170  make_degree_map(graph));
171  else
172  ::boost::cuthill_mckee_ordering(graph,
173  inv_perm.begin(),
174  get(::boost::vertex_color, graph),
175  make_degree_map(graph));
176 
177  for (boosttypes::size_type c = 0; c != inv_perm.size(); ++c)
178  new_dof_indices[index_map[inv_perm[c]]] = c;
179 
180  Assert(std::find(new_dof_indices.begin(),
181  new_dof_indices.end(),
182  numbers::invalid_dof_index) == new_dof_indices.end(),
183  ExcInternalError());
184  }
185 
186 
187 
188  template <typename DoFHandlerType>
189  void
190  king_ordering(DoFHandlerType &dof_handler,
191  const bool reversed_numbering,
192  const bool use_constraints)
193  {
194  std::vector<types::global_dof_index> renumbering(
195  dof_handler.n_dofs(), numbers::invalid_dof_index);
196  compute_king_ordering(renumbering,
197  dof_handler,
198  reversed_numbering,
199  use_constraints);
200 
201  // actually perform renumbering;
202  // this is dimension specific and
203  // thus needs an own function
204  dof_handler.renumber_dofs(renumbering);
205  }
206 
207 
208  template <typename DoFHandlerType>
209  void
210  compute_king_ordering(std::vector<types::global_dof_index> &new_dof_indices,
211  const DoFHandlerType & dof_handler,
212  const bool reversed_numbering,
213  const bool use_constraints)
214  {
215  boosttypes::Graph graph(dof_handler.n_dofs());
216  boosttypes::property_map<boosttypes::Graph,
217  boosttypes::vertex_degree_t>::type graph_degree;
218 
219  internal::create_graph(dof_handler, use_constraints, graph, graph_degree);
220 
221  boosttypes::property_map<boosttypes::Graph,
222  boosttypes::vertex_index_t>::type index_map =
223  get(::boost::vertex_index, graph);
224 
225 
226  std::vector<boosttypes::Vertex> inv_perm(num_vertices(graph));
227 
228  if (reversed_numbering == false)
229  ::boost::king_ordering(graph, inv_perm.rbegin());
230  else
231  ::boost::king_ordering(graph, inv_perm.begin());
232 
233  for (boosttypes::size_type c = 0; c != inv_perm.size(); ++c)
234  new_dof_indices[index_map[inv_perm[c]]] = c;
235 
236  Assert(std::find(new_dof_indices.begin(),
237  new_dof_indices.end(),
238  numbers::invalid_dof_index) == new_dof_indices.end(),
239  ExcInternalError());
240  }
241 
242 
243 
244  template <typename DoFHandlerType>
245  void
246  minimum_degree(DoFHandlerType &dof_handler,
247  const bool reversed_numbering,
248  const bool use_constraints)
249  {
250  std::vector<types::global_dof_index> renumbering(
251  dof_handler.n_dofs(), numbers::invalid_dof_index);
252  compute_minimum_degree(renumbering,
253  dof_handler,
254  reversed_numbering,
255  use_constraints);
256 
257  // actually perform renumbering;
258  // this is dimension specific and
259  // thus needs an own function
260  dof_handler.renumber_dofs(renumbering);
261  }
262 
263 
264  template <typename DoFHandlerType>
265  void
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)
271  {
272  (void)use_constraints;
273  Assert(use_constraints == false, ExcNotImplemented());
274 
275  // the following code is pretty
276  // much a verbatim copy of the
277  // sample code for the
278  // minimum_degree_ordering manual
279  // page from the BOOST Graph
280  // Library
281  using namespace ::boost;
282 
283  int delta = 0;
284 
285  // must be BGL directed graph now
286  using Graph = adjacency_list<vecS, vecS, directedS>;
287 
288  int n = dof_handler.n_dofs();
289 
290  Graph G(n);
291 
292  std::vector<::types::global_dof_index> dofs_on_this_cell;
293 
294  typename DoFHandlerType::active_cell_iterator cell = dof_handler
295  .begin_active(),
296  endc = dof_handler.end();
297 
298  for (; cell != endc; ++cell)
299  {
300  const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
301 
302  dofs_on_this_cell.resize(dofs_per_cell);
303 
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])
308  {
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);
311  }
312  }
313 
314 
315  using Vector = std::vector<int>;
316 
317 
318  Vector inverse_perm(n, 0);
319 
320  Vector perm(n, 0);
321 
322 
323  Vector supernode_sizes(n, 1);
324  // init has to be 1
325 
326  ::boost::property_map<Graph, vertex_index_t>::type id =
327  get(vertex_index, G);
328 
329 
330  Vector degree(n, 0);
331 
332 
333  minimum_degree_ordering(
334  G,
335  make_iterator_property_map(degree.begin(), id, degree[0]),
336  inverse_perm.data(),
337  perm.data(),
338  make_iterator_property_map(supernode_sizes.begin(),
339  id,
340  supernode_sizes[0]),
341  delta,
342  id);
343 
344 
345  for (int i = 0; i < n; ++i)
346  {
347  Assert(std::find(perm.begin(), perm.end(), i) != perm.end(),
348  ExcInternalError());
349  Assert(std::find(inverse_perm.begin(), inverse_perm.end(), i) !=
350  inverse_perm.end(),
351  ExcInternalError());
352  Assert(inverse_perm[perm[i]] == i, ExcInternalError());
353  }
354 
355  if (reversed_numbering == true)
356  std::copy(perm.begin(), perm.end(), new_dof_indices.begin());
357  else
358  std::copy(inverse_perm.begin(),
359  inverse_perm.end(),
360  new_dof_indices.begin());
361  }
362 
363  } // namespace boost
364 
365 
366 
367  template <typename DoFHandlerType>
368  void
369  Cuthill_McKee(DoFHandlerType & dof_handler,
370  const bool reversed_numbering,
371  const bool use_constraints,
372  const std::vector<types::global_dof_index> &starting_indices)
373  {
374  std::vector<types::global_dof_index> renumbering(
375  dof_handler.locally_owned_dofs().n_elements(),
377  compute_Cuthill_McKee(renumbering,
378  dof_handler,
379  reversed_numbering,
380  use_constraints,
381  starting_indices);
382 
383  // actually perform renumbering;
384  // this is dimension specific and
385  // thus needs an own function
386  dof_handler.renumber_dofs(renumbering);
387  }
388 
389 
390 
391  template <typename DoFHandlerType>
392  void
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)
399  {
400  // see if there is anything to do at all or whether we can skip the work on
401  // this processor
402  if (dof_handler.locally_owned_dofs().n_elements() == 0)
403  {
404  Assert(new_indices.size() == 0, ExcInternalError());
405  return;
406  }
407 
408  // make the connection graph
409  //
410  // note that if constraints are not requested, then the 'constraints'
411  // object will be empty and using it has no effect
412  IndexSet locally_relevant_dofs;
413  DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
414 
415  AffineConstraints<double> constraints;
416  if (use_constraints)
417  {
418  constraints.reinit(locally_relevant_dofs);
419  DoFTools::make_hanging_node_constraints(dof_handler, constraints);
420  }
421  constraints.close();
422 
423  const IndexSet &locally_owned_dofs = dof_handler.locally_owned_dofs();
424 
425  // see if we can get away with the sequential algorithm
426  if (locally_owned_dofs.n_elements() == locally_owned_dofs.size())
427  {
428  AssertDimension(new_indices.size(), dof_handler.n_dofs());
429 
430  DynamicSparsityPattern dsp(dof_handler.n_dofs(), dof_handler.n_dofs());
431  DoFTools::make_sparsity_pattern(dof_handler, dsp, constraints);
432 
434  new_indices,
435  starting_indices);
436  if (reversed_numbering)
437  new_indices = Utilities::reverse_permutation(new_indices);
438  }
439  else
440  {
441  // we are in the parallel case where we need to work in the
442  // local index space, i.e., the locally owned part of the
443  // sparsity pattern.
444  //
445  // first figure out whether the user only gave us starting
446  // indices that are locally owned, or that are only locally
447  // relevant. in the process, also check that all indices
448  // really belong to at least the locally relevant ones
449  IndexSet locally_active_dofs;
450  DoFTools::extract_locally_active_dofs(dof_handler, locally_active_dofs);
451 
452  bool needs_locally_active = false;
453  for (const auto starting_index : starting_indices)
454  {
455  if ((needs_locally_active ==
456  /* previously already set to */ true) ||
457  (locally_owned_dofs.is_element(starting_index) == false))
458  {
459  Assert(
460  locally_active_dofs.is_element(starting_index),
461  ExcMessage(
462  "You specified global degree of freedom " +
463  std::to_string(starting_index) +
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;
468  }
469  }
470 
471  const IndexSet index_set_to_use =
472  (needs_locally_active ? locally_active_dofs : locally_owned_dofs);
473 
474  // then create first the global sparsity pattern, and then the local
475  // sparsity pattern from the global one by transferring its indices to
476  // processor-local (locally owned or locally active) index space
477  DynamicSparsityPattern dsp(dof_handler.n_dofs(),
478  dof_handler.n_dofs(),
479  index_set_to_use);
480  DoFTools::make_sparsity_pattern(dof_handler, dsp, constraints);
481 
482  DynamicSparsityPattern local_sparsity(index_set_to_use.n_elements(),
483  index_set_to_use.n_elements());
484  std::vector<types::global_dof_index> row_entries;
485  for (unsigned int i = 0; i < index_set_to_use.n_elements(); ++i)
486  {
487  const types::global_dof_index row =
488  index_set_to_use.nth_index_in_set(i);
489  const unsigned int row_length = dsp.row_length(row);
490  row_entries.clear();
491  for (unsigned int j = 0; j < row_length; ++j)
492  {
493  const unsigned int col = dsp.column_number(row, j);
494  if (col != row && index_set_to_use.is_element(col))
495  row_entries.push_back(index_set_to_use.index_within_set(col));
496  }
497  local_sparsity.add_entries(i,
498  row_entries.begin(),
499  row_entries.end(),
500  true);
501  }
502 
503  // translate starting indices from global to local indices
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] =
508  index_set_to_use.index_within_set(starting_indices[i]);
509 
510  // then do the renumbering on the locally owned portion
511  AssertDimension(new_indices.size(), locally_owned_dofs.n_elements());
512  std::vector<types::global_dof_index> my_new_indices(
513  index_set_to_use.n_elements());
515  my_new_indices,
516  local_starting_indices);
517  if (reversed_numbering)
518  my_new_indices = Utilities::reverse_permutation(my_new_indices);
519 
520  // now that we have a re-enumeration of all DoFs, we need to throw
521  // out the ones that are not locally owned in case we have worked
522  // with the locally active ones. that's because the renumbering
523  // functions only want new indices for the locally owned DoFs (other
524  // processors are responsible for renumbering the ones that are
525  // on cell interfaces)
526  if (needs_locally_active == true)
527  {
528  // first step: figure out which DoF indices to eliminate
529  IndexSet active_but_not_owned_dofs = locally_active_dofs;
530  active_but_not_owned_dofs.subtract_set(locally_owned_dofs);
531 
532  std::set<types::global_dof_index> erase_these_indices;
533  for (const auto p : active_but_not_owned_dofs)
534  {
535  const auto index = index_set_to_use.index_within_set(p);
536  Assert(index < index_set_to_use.n_elements(),
537  ExcInternalError());
538  erase_these_indices.insert(my_new_indices[index]);
539  my_new_indices[index] = numbers::invalid_dof_index;
540  }
541  Assert(erase_these_indices.size() ==
542  active_but_not_owned_dofs.n_elements(),
543  ExcInternalError());
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(),
549  ExcInternalError());
550 
551  // then compute a renumbering of the remaining ones
552  std::vector<types::global_dof_index> translate_indices(
553  my_new_indices.size());
554  {
555  std::set<types::global_dof_index>::const_iterator
556  next_erased_index = erase_these_indices.begin();
557  types::global_dof_index next_new_index = 0;
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))
561  {
562  translate_indices[i] = numbers::invalid_dof_index;
563  ++next_erased_index;
564  }
565  else
566  {
567  translate_indices[i] = next_new_index;
568  ++next_new_index;
569  }
570  Assert(next_new_index == locally_owned_dofs.n_elements(),
571  ExcInternalError());
572  }
573 
574  // and then do the renumbering of the result of the
575  // Cuthill-McKee algorithm above, right into the output array
576  new_indices.clear();
577  new_indices.reserve(locally_owned_dofs.n_elements());
578  for (const auto &p : my_new_indices)
580  {
581  Assert(translate_indices[p] != numbers::invalid_dof_index,
582  ExcInternalError());
583  new_indices.push_back(translate_indices[p]);
584  }
585  Assert(new_indices.size() == locally_owned_dofs.n_elements(),
586  ExcInternalError());
587  }
588  else
589  new_indices = std::move(my_new_indices);
590 
591  // convert indices back to global index space. in both of the branches
592  // above, we ended up with new_indices only containing the local
593  // indices of the locally-owned DoFs. so that's where we get the
594  // indices
595  for (types::global_dof_index &new_index : new_indices)
596  new_index = locally_owned_dofs.nth_index_in_set(new_index);
597  }
598  }
599 
600 
601 
602  template <typename DoFHandlerType>
603  void
604  Cuthill_McKee(DoFHandlerType & dof_handler,
605  const unsigned int level,
606  const bool reversed_numbering,
607  const std::vector<types::global_dof_index> &starting_indices)
608  {
609  Assert(dof_handler.n_dofs(level) != numbers::invalid_dof_index,
611 
612  // make the connection graph
613  DynamicSparsityPattern dsp(dof_handler.n_dofs(level),
614  dof_handler.n_dofs(level));
615  MGTools::make_sparsity_pattern(dof_handler, dsp, level);
616 
617  std::vector<types::global_dof_index> new_indices(dsp.n_rows());
618  SparsityTools::reorder_Cuthill_McKee(dsp, new_indices, starting_indices);
619 
620  if (reversed_numbering)
621  new_indices = Utilities::reverse_permutation(new_indices);
622 
623  // actually perform renumbering;
624  // this is dimension specific and
625  // thus needs an own function
626  dof_handler.renumber_dofs(level, new_indices);
627  }
628 
629 
630 
631  template <typename DoFHandlerType>
632  void
633  component_wise(DoFHandlerType & dof_handler,
634  const std::vector<unsigned int> &component_order_arg)
635  {
636  std::vector<types::global_dof_index> renumbering(
637  dof_handler.n_locally_owned_dofs(), numbers::invalid_dof_index);
638 
639  typename DoFHandlerType::active_cell_iterator start =
640  dof_handler.begin_active();
641  const typename DoFHandlerType::level_cell_iterator end = dof_handler.end();
642 
643  const types::global_dof_index result =
644  compute_component_wise<DoFHandlerType::dimension,
645  DoFHandlerType::space_dimension>(
646  renumbering, start, end, component_order_arg, false);
647  if (result == 0)
648  return;
649 
650  // verify that the last numbered
651  // degree of freedom is either
652  // equal to the number of degrees
653  // of freedom in total (the
654  // sequential case) or in the
655  // distributed case at least
656  // makes sense
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())),
660  ExcInternalError());
661 
662  dof_handler.renumber_dofs(renumbering);
663  }
664 
665 
666 
667  template <typename DoFHandlerType>
668  void
669  component_wise(DoFHandlerType & dof_handler,
670  const unsigned int level,
671  const std::vector<unsigned int> &component_order_arg)
672  {
673  Assert(dof_handler.n_dofs(level) != numbers::invalid_dof_index,
675 
676  std::vector<types::global_dof_index> renumbering(
677  dof_handler.locally_owned_mg_dofs(level).n_elements(),
679 
680  typename DoFHandlerType::level_cell_iterator start =
681  dof_handler.begin(level);
682  typename DoFHandlerType::level_cell_iterator end = dof_handler.end(level);
683 
684  const types::global_dof_index result =
685  compute_component_wise<DoFHandlerType::dimension,
686  DoFHandlerType::space_dimension>(
687  renumbering, start, end, component_order_arg, true);
688 
689  if (result == 0)
690  return;
691 
692  Assert(result == dof_handler.n_dofs(level), ExcInternalError());
693 
694  if (renumbering.size() != 0)
695  dof_handler.renumber_dofs(level, renumbering);
696  }
697 
698 
699 
700  template <int dim, int spacedim, typename CellIterator>
702  compute_component_wise(std::vector<types::global_dof_index> &new_indices,
703  const CellIterator & start,
704  const typename identity<CellIterator>::type &end,
705  const std::vector<unsigned int> &component_order_arg,
706  const bool is_level_operation)
707  {
708  const hp::FECollection<dim, spacedim> fe_collection(
709  start->get_dof_handler().get_fe_collection());
710 
711  // do nothing if the FE has only
712  // one component
713  if (fe_collection.n_components() == 1)
714  {
715  new_indices.resize(0);
716  return 0;
717  }
718 
719  // Get a reference to the set of dofs. Note that we assume that all cells
720  // are assumed to be on the same level, otherwise the operation doesn't make
721  // much sense (we will assert this below).
722  const IndexSet &locally_owned_dofs =
723  is_level_operation ?
724  start->get_dof_handler().locally_owned_mg_dofs(start->level()) :
725  start->get_dof_handler().locally_owned_dofs();
726 
727  // Copy last argument into a
728  // writable vector.
729  std::vector<unsigned int> component_order(component_order_arg);
730  // If the last argument was an
731  // empty vector, set up things to
732  // store components in the order
733  // found in the system.
734  if (component_order.size() == 0)
735  for (unsigned int i = 0; i < fe_collection.n_components(); ++i)
736  component_order.push_back(i);
737 
738  Assert(component_order.size() == fe_collection.n_components(),
739  ExcDimensionMismatch(component_order.size(),
740  fe_collection.n_components()));
741 
742  for (const unsigned int component : component_order)
743  {
744  (void)component;
745  AssertIndexRange(component, fe_collection.n_components());
746  }
747 
748  // vector to hold the dof indices on
749  // the cell we visit at a time
750  std::vector<types::global_dof_index> local_dof_indices;
751 
752  // prebuilt list to which component
753  // a given dof on a cell
754  // should go. note that we get into
755  // trouble here if the shape
756  // function is not primitive, since
757  // then there is no single vector
758  // component to which it
759  // belongs. in this case, assign it
760  // to the first vector component to
761  // which it belongs
762  std::vector<std::vector<unsigned int>> component_list(fe_collection.size());
763  for (unsigned int f = 0; f < fe_collection.size(); ++f)
764  {
765  const FiniteElement<dim, spacedim> &fe = fe_collection[f];
766  const unsigned int dofs_per_cell = fe.dofs_per_cell;
767  component_list[f].resize(dofs_per_cell);
768  for (unsigned int i = 0; i < dofs_per_cell; ++i)
769  if (fe.is_primitive(i))
770  component_list[f][i] =
771  component_order[fe.system_to_component_index(i).first];
772  else
773  {
774  const unsigned int comp =
776 
777  // then associate this degree
778  // of freedom with this
779  // component
780  component_list[f][i] = component_order[comp];
781  }
782  }
783 
784  // set up a map where for each
785  // component the respective degrees
786  // of freedom are collected.
787  //
788  // note that this map is sorted by
789  // component but that within each
790  // component it is NOT sorted by
791  // dof index. note also that some
792  // dof indices are entered
793  // multiply, so we will have to
794  // take care of that
795  std::vector<std::vector<types::global_dof_index>> component_to_dof_map(
796  fe_collection.n_components());
797  for (CellIterator cell = start; cell != end; ++cell)
798  {
799  if (is_level_operation)
800  {
801  // we are dealing with mg dofs, skip foreign level cells:
802  if (!cell->is_locally_owned_on_level())
803  continue;
804  }
805  else
806  {
807  // we are dealing with active dofs, skip the loop if not locally
808  // owned:
809  if (!cell->is_locally_owned())
810  continue;
811  }
812 
813  if (is_level_operation)
814  Assert(
815  cell->level() == start->level(),
816  ExcMessage(
817  "Multigrid renumbering in compute_component_wise() needs to be applied to a single level!"));
818 
819  // on each cell: get dof indices
820  // and insert them into the global
821  // list using their component
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);
827 
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]);
832  }
833 
834  // now we've got all indices sorted
835  // into buckets labeled by their
836  // target component number. we've
837  // only got to traverse this list
838  // and assign the new indices
839  //
840  // however, we first want to sort
841  // the indices entered into the
842  // buckets to preserve the order
843  // within each component and during
844  // this also remove duplicate
845  // entries
846  //
847  // note that we no longer have to
848  // care about non-primitive shape
849  // functions since the buckets
850  // corresponding to the second and
851  // following vector components of a
852  // non-primitive FE will simply be
853  // empty, everything being shoved
854  // into the first one. The same
855  // holds if several components were
856  // joined into a single target.
857  for (unsigned int component = 0; component < fe_collection.n_components();
858  ++component)
859  {
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());
866  }
867 
868  // calculate the number of locally owned
869  // DoFs per bucket
870  const unsigned int n_buckets = fe_collection.n_components();
871  std::vector<types::global_dof_index> shifts(n_buckets);
872 
874  (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
875  &start->get_dof_handler().get_triangulation())))
876  {
877 #ifdef DEAL_II_WITH_MPI
878  std::vector<types::global_dof_index> local_dof_count(n_buckets);
879 
880  for (unsigned int c = 0; c < n_buckets; ++c)
881  local_dof_count[c] = component_to_dof_map[c].size();
882 
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(),
886  n_buckets,
888  MPI_SUM,
889  tria->get_communicator());
890  AssertThrowMPI(ierr);
891 
892  std::vector<types::global_dof_index> global_dof_count(n_buckets);
893  Utilities::MPI::sum(local_dof_count,
894  tria->get_communicator(),
895  global_dof_count);
896 
897  // calculate shifts
898  types::global_dof_index cumulated = 0;
899  for (unsigned int c = 0; c < n_buckets; ++c)
900  {
901  shifts[c] = prefix_dof_count[c] + cumulated;
902  cumulated += global_dof_count[c];
903  }
904 #else
905  (void)tria;
906  Assert(false, ExcInternalError());
907 #endif
908  }
909  else
910  {
911  shifts[0] = 0;
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();
914  }
915 
916 
917 
918  // now concatenate all the
919  // components in the order the user
920  // desired to see
921  types::global_dof_index next_free_index = 0;
922  for (unsigned int component = 0; component < fe_collection.n_components();
923  ++component)
924  {
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();
928 
929  next_free_index = shifts[component];
930 
931  for (typename std::vector<types::global_dof_index>::const_iterator
932  dof_index = begin_of_component;
933  dof_index != end_of_component;
934  ++dof_index)
935  {
936  Assert(locally_owned_dofs.index_within_set(*dof_index) <
937  new_indices.size(),
938  ExcInternalError());
939  new_indices[locally_owned_dofs.index_within_set(*dof_index)] =
940  next_free_index++;
941  }
942  }
943 
944  return next_free_index;
945  }
946 
947 
948 
949  template <int dim, int spacedim>
950  void
952  {
953  std::vector<types::global_dof_index> renumbering(
955 
957  dof_handler.begin_active();
959  dof_handler.end();
960 
962  dim,
963  spacedim,
966  start,
967  end,
968  false);
969  if (result == 0)
970  return;
971 
972  // verify that the last numbered
973  // degree of freedom is either
974  // equal to the number of degrees
975  // of freedom in total (the
976  // sequential case) or in the
977  // distributed case at least
978  // makes sense
979  Assert((result == dof_handler.n_locally_owned_dofs()) ||
980  ((dof_handler.n_locally_owned_dofs() < dof_handler.n_dofs()) &&
981  (result <= dof_handler.n_dofs())),
982  ExcInternalError());
983 
984  dof_handler.renumber_dofs(renumbering);
985  }
986 
987 
988 
989  template <int dim, int spacedim>
990  void
992  {
993  std::vector<types::global_dof_index> renumbering(
994  dof_handler.n_dofs(), numbers::invalid_dof_index);
995 
997  dof_handler.begin_active();
999  dof_handler.end();
1000 
1002  dim,
1003  spacedim,
1006  start,
1007  end,
1008  false);
1009 
1010  if (result == 0)
1011  return;
1012 
1013  Assert(result == dof_handler.n_dofs(), ExcInternalError());
1014 
1015  dof_handler.renumber_dofs(renumbering);
1016  }
1017 
1018 
1019 
1020  template <int dim, int spacedim>
1021  void
1022  block_wise(DoFHandler<dim, spacedim> &dof_handler, const unsigned int level)
1023  {
1026 
1027  std::vector<types::global_dof_index> renumbering(
1028  dof_handler.locally_owned_mg_dofs(level).n_elements(),
1030 
1032  dof_handler.begin(level);
1034  dof_handler.end(level);
1035 
1037  dim,
1038  spacedim,
1040  typename DoFHandler<dim, spacedim>::level_cell_iterator>(renumbering,
1041  start,
1042  end,
1043  true);
1044 
1045  if (result == 0)
1046  return;
1047 
1048  Assert(result == dof_handler.n_dofs(level), ExcInternalError());
1049 
1050  if (renumbering.size() != 0)
1051  dof_handler.renumber_dofs(level, renumbering);
1052  }
1053 
1054 
1055 
1056  template <int dim, int spacedim, class ITERATOR, class ENDITERATOR>
1058  compute_block_wise(std::vector<types::global_dof_index> &new_indices,
1059  const ITERATOR & start,
1060  const ENDITERATOR & end,
1061  const bool is_level_operation)
1062  {
1063  const hp::FECollection<dim, spacedim> fe_collection(
1064  start->get_dof_handler().get_fe_collection());
1065 
1066  // do nothing if the FE has only
1067  // one component
1068  if (fe_collection.n_blocks() == 1)
1069  {
1070  new_indices.resize(0);
1071  return 0;
1072  }
1073 
1074  // Get a reference to the set of dofs. Note that we assume that all cells
1075  // are assumed to be on the same level, otherwise the operation doesn't make
1076  // much sense (we will assert this below).
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();
1081 
1082  // vector to hold the dof indices on
1083  // the cell we visit at a time
1084  std::vector<types::global_dof_index> local_dof_indices;
1085 
1086  // prebuilt list to which block
1087  // a given dof on a cell
1088  // should go.
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)
1092  {
1093  const FiniteElement<dim, spacedim> &fe = fe_collection[f];
1094  block_list[f].resize(fe.dofs_per_cell);
1095  for (unsigned int i = 0; i < fe.dofs_per_cell; ++i)
1096  block_list[f][i] = fe.system_to_block_index(i).first;
1097  }
1098 
1099  // set up a map where for each
1100  // block the respective degrees
1101  // of freedom are collected.
1102  //
1103  // note that this map is sorted by
1104  // block but that within each
1105  // block it is NOT sorted by
1106  // dof index. note also that some
1107  // dof indices are entered
1108  // multiply, so we will have to
1109  // take care of that
1110  std::vector<std::vector<types::global_dof_index>> block_to_dof_map(
1111  fe_collection.n_blocks());
1112  for (ITERATOR cell = start; cell != end; ++cell)
1113  {
1114  if (is_level_operation)
1115  {
1116  // we are dealing with mg dofs, skip foreign level cells:
1117  if (!cell->is_locally_owned_on_level())
1118  continue;
1119  }
1120  else
1121  {
1122  // we are dealing with active dofs, skip the loop if not locally
1123  // owned:
1124  if (!cell->is_locally_owned())
1125  continue;
1126  }
1127 
1128  if (is_level_operation)
1129  Assert(
1130  cell->level() == start->level(),
1131  ExcMessage(
1132  "Multigrid renumbering in compute_block_wise() needs to be applied to a single level!"));
1133 
1134  // on each cell: get dof indices
1135  // and insert them into the global
1136  // list using their component
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);
1142 
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]);
1147  }
1148 
1149  // now we've got all indices sorted
1150  // into buckets labeled by their
1151  // target block number. we've
1152  // only got to traverse this list
1153  // and assign the new indices
1154  //
1155  // however, we first want to sort
1156  // the indices entered into the
1157  // buckets to preserve the order
1158  // within each component and during
1159  // this also remove duplicate
1160  // entries
1161  for (unsigned int block = 0; block < fe_collection.n_blocks(); ++block)
1162  {
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());
1169  }
1170 
1171  // calculate the number of locally owned
1172  // DoFs per bucket
1173  const unsigned int n_buckets = fe_collection.n_blocks();
1174  std::vector<types::global_dof_index> shifts(n_buckets);
1175 
1177  (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1178  &start->get_dof_handler().get_triangulation())))
1179  {
1180 #ifdef DEAL_II_WITH_MPI
1181  std::vector<types::global_dof_index> local_dof_count(n_buckets);
1182 
1183  for (unsigned int c = 0; c < n_buckets; ++c)
1184  local_dof_count[c] = block_to_dof_map[c].size();
1185 
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(),
1189  n_buckets,
1191  MPI_SUM,
1192  tria->get_communicator());
1193  AssertThrowMPI(ierr);
1194 
1195  std::vector<types::global_dof_index> global_dof_count(n_buckets);
1196  Utilities::MPI::sum(local_dof_count,
1197  tria->get_communicator(),
1198  global_dof_count);
1199 
1200  // calculate shifts
1201  types::global_dof_index cumulated = 0;
1202  for (unsigned int c = 0; c < n_buckets; ++c)
1203  {
1204  shifts[c] = prefix_dof_count[c] + cumulated;
1205  cumulated += global_dof_count[c];
1206  }
1207 #else
1208  (void)tria;
1209  Assert(false, ExcInternalError());
1210 #endif
1211  }
1212  else
1213  {
1214  shifts[0] = 0;
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();
1217  }
1218 
1219 
1220 
1221  // now concatenate all the
1222  // components in the order the user
1223  // desired to see
1224  types::global_dof_index next_free_index = 0;
1225  for (unsigned int block = 0; block < fe_collection.n_blocks(); ++block)
1226  {
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();
1230 
1231  next_free_index = shifts[block];
1232 
1233  for (typename std::vector<types::global_dof_index>::const_iterator
1234  dof_index = begin_of_component;
1235  dof_index != end_of_component;
1236  ++dof_index)
1237  {
1238  Assert(locally_owned_dofs.index_within_set(*dof_index) <
1239  new_indices.size(),
1240  ExcInternalError());
1241  new_indices[locally_owned_dofs.index_within_set(*dof_index)] =
1242  next_free_index++;
1243  }
1244  }
1245 
1246  return next_free_index;
1247  }
1248 
1249 
1250 
1251  namespace
1252  {
1253  // Helper function for DoFRenumbering::hierarchical(). This function
1254  // recurses into the given cell or, if that should be an active (terminal)
1255  // cell, renumbers DoF indices on it. The function starts renumbering with
1256  // 'next_free_dof_index' and returns the first still unused DoF index at the
1257  // end of its operation.
1258  template <int dim, class CellIteratorType>
1260  compute_hierarchical_recursive(
1261  const types::global_dof_index next_free_dof_offset,
1262  const types::global_dof_index my_starting_index,
1263  const CellIteratorType & cell,
1264  const IndexSet & locally_owned_dof_indices,
1265  std::vector<types::global_dof_index> &new_indices)
1266  {
1267  types::global_dof_index current_next_free_dof_offset =
1268  next_free_dof_offset;
1269 
1270  if (cell->has_children())
1271  {
1272  // recursion
1273  for (unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell;
1274  ++c)
1275  current_next_free_dof_offset =
1276  compute_hierarchical_recursive<dim>(current_next_free_dof_offset,
1277  my_starting_index,
1278  cell->child(c),
1279  locally_owned_dof_indices,
1280  new_indices);
1281  }
1282  else
1283  {
1284  // this is a terminal cell. we need to renumber its DoF indices. there
1285  // are now three cases to decide:
1286  // - this is a sequential triangulation: we can just go ahead and
1287  // number
1288  // the DoFs in the order in which we encounter cells. in this case,
1289  // all cells are actually locally owned
1290  // - if this is a parallel::distributed::Triangulation, then we only
1291  // need to work on the locally owned cells since they contain
1292  // all locally owned DoFs.
1293  // - if this is a parallel::shared::Triangulation, then the same
1294  // applies
1295  //
1296  // in all cases, each processor starts new indices so that we get
1297  // a consecutive numbering on each processor, and disjoint ownership
1298  // of the global range of DoF indices
1299  if (cell->is_locally_owned())
1300  {
1301  // first get the existing DoF indices
1302  const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
1303  std::vector<types::global_dof_index> local_dof_indices(
1304  dofs_per_cell);
1305  cell->get_dof_indices(local_dof_indices);
1306 
1307  // then loop over the existing DoF indices on this cell
1308  // and see whether it has already been re-numbered (it
1309  // may have been on a face or vertex to a neighboring
1310  // cell that we have encountered before). if not,
1311  // give it a new number and store that number in the
1312  // output array (new_indices)
1313  //
1314  // if this is a parallel triangulation and a DoF index is
1315  // not locally owned, then don't touch it. since
1316  // we don't actually *change* DoF indices (just record new
1317  // numbers in an array), we don't need to worry about
1318  // the decision whether a DoF is locally owned or not changing
1319  // as we progress in renumbering DoFs -- all adjacent cells
1320  // will always agree that a DoF is locally owned or not.
1321  // that said, the first cell to encounter a locally owned DoF
1322  // gets to number it, so the order in which we traverse cells
1323  // matters
1324  for (unsigned int i = 0; i < dofs_per_cell; ++i)
1325  if (locally_owned_dof_indices.is_element(local_dof_indices[i]))
1326  {
1327  // this is a locally owned DoF, assign new number if not
1328  // assigned a number yet
1329  const unsigned int idx =
1330  locally_owned_dof_indices.index_within_set(
1331  local_dof_indices[i]);
1332  if (new_indices[idx] == numbers::invalid_dof_index)
1333  {
1334  new_indices[idx] =
1335  my_starting_index + current_next_free_dof_offset;
1336  ++current_next_free_dof_offset;
1337  }
1338  }
1339  }
1340  }
1341 
1342  return current_next_free_dof_offset;
1343  }
1344  } // namespace
1345 
1346 
1347 
1348  template <typename DoFHandlerType>
1349  void
1350  hierarchical(DoFHandlerType &dof_handler)
1351  {
1352  const int dim = DoFHandlerType::dimension;
1353  const int spacedim = DoFHandlerType::space_dimension;
1354 
1355  std::vector<types::global_dof_index> renumbering(
1356  dof_handler.n_locally_owned_dofs(), numbers::invalid_dof_index);
1357 
1358  types::global_dof_index next_free_dof_offset = 0;
1359  const IndexSet locally_owned = dof_handler.locally_owned_dofs();
1360 
1361  // in the function we call recursively, we want to number DoFs so
1362  // that global cell zero starts with DoF zero, regardless of how
1363  // DoFs were previously numbered. to this end, we need to figure
1364  // out which DoF index the current processor should start with.
1365  //
1366  // if this is a sequential triangulation, then obviously the starting
1367  // index is zero. otherwise, make sure we get contiguous, successive
1368  // ranges on each processor. note that since the number of locally owned
1369  // DoFs is an invariant under renumbering, we can easily compute this
1370  // starting index by just accumulating over the number of locally owned
1371  // DoFs for all previous processes
1372  types::global_dof_index my_starting_index = 0;
1373 
1375  dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1376  &dof_handler.get_triangulation()))
1377  {
1378 #ifdef DEAL_II_WITH_MPI
1379  types::global_dof_index local_size =
1380  dof_handler.locally_owned_dofs().n_elements();
1381  MPI_Exscan(&local_size,
1382  &my_starting_index,
1383  1,
1385  MPI_SUM,
1386  tria->get_communicator());
1387 #endif
1388  }
1389 
1392  *>(&dof_handler.get_triangulation()))
1393  {
1394 #ifdef DEAL_II_WITH_P4EST
1395  // this is a distributed Triangulation. we need to traverse the coarse
1396  // cells in the order p4est does to match the z-order actually used
1397  // by p4est. this requires using the renumbering of coarse cells
1398  // we do before we hand things off to p4est
1399  for (unsigned int c = 0; c < tria->n_cells(0); ++c)
1400  {
1401  const unsigned int coarse_cell_index =
1402  tria->get_p4est_tree_to_coarse_cell_permutation()[c];
1403 
1404  const typename DoFHandlerType::level_cell_iterator this_cell(
1405  tria, 0, coarse_cell_index, &dof_handler);
1406 
1407  next_free_dof_offset =
1408  compute_hierarchical_recursive<dim>(next_free_dof_offset,
1409  my_starting_index,
1410  this_cell,
1411  locally_owned,
1412  renumbering);
1413  }
1414 #else
1415  Assert(false, ExcNotImplemented());
1416 #endif
1417  }
1418  else
1419  {
1420  // this is not a distributed Triangulation, so we can traverse coarse
1421  // cells in the normal order
1422  for (typename DoFHandlerType::cell_iterator cell = dof_handler.begin(0);
1423  cell != dof_handler.end(0);
1424  ++cell)
1425  next_free_dof_offset =
1426  compute_hierarchical_recursive<dim>(next_free_dof_offset,
1427  my_starting_index,
1428  cell,
1429  locally_owned,
1430  renumbering);
1431  }
1432 
1433  // verify that the last numbered degree of freedom is either
1434  // equal to the number of degrees of freedom in total (the
1435  // sequential case) or in the distributed case at least
1436  // makes sense
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())),
1440  ExcInternalError());
1441 
1442  // make sure that all local DoFs got new numbers assigned
1443  Assert(std::find(renumbering.begin(),
1444  renumbering.end(),
1445  numbers::invalid_dof_index) == renumbering.end(),
1446  ExcInternalError());
1447 
1448  dof_handler.renumber_dofs(renumbering);
1449  }
1450 
1451 
1452 
1453  template <typename DoFHandlerType>
1454  void
1455  sort_selected_dofs_back(DoFHandlerType & dof_handler,
1456  const std::vector<bool> &selected_dofs)
1457  {
1458  std::vector<types::global_dof_index> renumbering(
1459  dof_handler.n_dofs(), numbers::invalid_dof_index);
1460  compute_sort_selected_dofs_back(renumbering, dof_handler, selected_dofs);
1461 
1462  dof_handler.renumber_dofs(renumbering);
1463  }
1464 
1465 
1466 
1467  template <typename DoFHandlerType>
1468  void
1469  sort_selected_dofs_back(DoFHandlerType & dof_handler,
1470  const std::vector<bool> &selected_dofs,
1471  const unsigned int level)
1472  {
1473  Assert(dof_handler.n_dofs(level) != numbers::invalid_dof_index,
1475 
1476  std::vector<types::global_dof_index> renumbering(
1477  dof_handler.n_dofs(level), numbers::invalid_dof_index);
1478  compute_sort_selected_dofs_back(renumbering,
1479  dof_handler,
1480  selected_dofs,
1481  level);
1482 
1483  dof_handler.renumber_dofs(level, renumbering);
1484  }
1485 
1486 
1487 
1488  template <typename DoFHandlerType>
1489  void
1491  std::vector<types::global_dof_index> &new_indices,
1492  const DoFHandlerType & dof_handler,
1493  const std::vector<bool> & selected_dofs)
1494  {
1495  const types::global_dof_index n_dofs = dof_handler.n_dofs();
1496  Assert(selected_dofs.size() == n_dofs,
1497  ExcDimensionMismatch(selected_dofs.size(), n_dofs));
1498 
1499  // re-sort the dofs according to
1500  // their selection state
1501  Assert(new_indices.size() == n_dofs,
1502  ExcDimensionMismatch(new_indices.size(), n_dofs));
1503 
1504  const types::global_dof_index n_selected_dofs =
1505  std::count(selected_dofs.begin(), selected_dofs.end(), false);
1506 
1507  types::global_dof_index next_unselected = 0;
1508  types::global_dof_index next_selected = n_selected_dofs;
1509  for (types::global_dof_index i = 0; i < n_dofs; ++i)
1510  if (selected_dofs[i] == false)
1511  {
1512  new_indices[i] = next_unselected;
1513  ++next_unselected;
1514  }
1515  else
1516  {
1517  new_indices[i] = next_selected;
1518  ++next_selected;
1519  }
1520  Assert(next_unselected == n_selected_dofs, ExcInternalError());
1521  Assert(next_selected == n_dofs, ExcInternalError());
1522  }
1523 
1524 
1525 
1526  template <typename DoFHandlerType>
1527  void
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)
1533  {
1534  Assert(dof_handler.n_dofs(level) != numbers::invalid_dof_index,
1536 
1537  const unsigned int n_dofs = dof_handler.n_dofs(level);
1538  Assert(selected_dofs.size() == n_dofs,
1539  ExcDimensionMismatch(selected_dofs.size(), n_dofs));
1540 
1541  // re-sort the dofs according to
1542  // their selection state
1543  Assert(new_indices.size() == n_dofs,
1544  ExcDimensionMismatch(new_indices.size(), n_dofs));
1545 
1546  const unsigned int n_selected_dofs =
1547  std::count(selected_dofs.begin(), selected_dofs.end(), false);
1548 
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)
1553  {
1554  new_indices[i] = next_unselected;
1555  ++next_unselected;
1556  }
1557  else
1558  {
1559  new_indices[i] = next_selected;
1560  ++next_selected;
1561  }
1562  Assert(next_unselected == n_selected_dofs, ExcInternalError());
1563  Assert(next_selected == n_dofs, ExcInternalError());
1564  }
1565 
1566 
1567 
1568  template <typename DoFHandlerType>
1569  void
1571  DoFHandlerType & dof,
1572  const std::vector<typename DoFHandlerType::active_cell_iterator> &cells)
1573  {
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());
1577  compute_cell_wise(renumbering, reverse, dof, cells);
1578 
1579  dof.renumber_dofs(renumbering);
1580  }
1581 
1582 
1583  template <typename DoFHandlerType>
1584  void
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>
1590  &cells)
1591  {
1592  if (const parallel::TriangulationBase<DoFHandlerType::dimension,
1593  DoFHandlerType::space_dimension> *p =
1594  dynamic_cast<const parallel::TriangulationBase<
1595  DoFHandlerType::dimension,
1596  DoFHandlerType::space_dimension> *>(&dof.get_triangulation()))
1597  {
1598  AssertDimension(cells.size(), p->n_locally_owned_active_cells());
1599  }
1600  else
1601  {
1602  AssertDimension(cells.size(), dof.get_triangulation().n_active_cells());
1603  }
1604 
1605  const auto n_owned_dofs = dof.n_locally_owned_dofs();
1606 
1607  // Actually, we compute the inverse of the reordering vector, called reverse
1608  // here. Later, its inverse is computed into new_indices, which is the
1609  // return argument.
1610 
1611  AssertDimension(new_indices.size(), n_owned_dofs);
1612  AssertDimension(reverse.size(), n_owned_dofs);
1613 
1614  // For continuous elements, we must make sure, that each dof is reordered
1615  // only once.
1616  std::vector<bool> already_sorted(n_owned_dofs, false);
1617  std::vector<types::global_dof_index> cell_dofs;
1618 
1619  const auto &owned_dofs = dof.locally_owned_dofs();
1620 
1621  unsigned int index = 0;
1622 
1623  for (const auto &cell : cells)
1624  {
1625  // Determine the number of dofs on this cell and reinit the
1626  // vector storing these numbers.
1627  const unsigned int n_cell_dofs = cell->get_fe().n_dofs_per_cell();
1628  cell_dofs.resize(n_cell_dofs);
1629 
1630  cell->get_active_or_mg_dof_indices(cell_dofs);
1631 
1632  // Sort here to make sure that degrees of freedom inside a single cell
1633  // are in the same order after renumbering.
1634  std::sort(cell_dofs.begin(), cell_dofs.end());
1635 
1636  for (const auto dof : cell_dofs)
1637  {
1638  const auto local_dof = owned_dofs.index_within_set(dof);
1639  if (local_dof != numbers::invalid_dof_index &&
1640  !already_sorted[local_dof])
1641  {
1642  already_sorted[local_dof] = true;
1643  reverse[index++] = local_dof;
1644  }
1645  }
1646  }
1647  Assert(index == n_owned_dofs,
1648  ExcMessage(
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?"));
1652 
1653  for (types::global_dof_index i = 0; i < reverse.size(); ++i)
1654  new_indices[reverse[i]] = owned_dofs.nth_index_in_set(i);
1655  }
1656 
1657 
1658 
1659  template <typename DoFHandlerType>
1660  void
1662  DoFHandlerType & dof,
1663  const unsigned int level,
1664  const typename std::vector<typename DoFHandlerType::level_cell_iterator>
1665  &cells)
1666  {
1667  Assert(dof.n_dofs(level) != numbers::invalid_dof_index,
1669 
1670  std::vector<types::global_dof_index> renumbering(dof.n_dofs(level));
1671  std::vector<types::global_dof_index> reverse(dof.n_dofs(level));
1672 
1673  compute_cell_wise(renumbering, reverse, dof, level, cells);
1674  dof.renumber_dofs(level, renumbering);
1675  }
1676 
1677 
1678 
1679  template <typename DoFHandlerType>
1680  void
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>
1687  &cells)
1688  {
1689  Assert(cells.size() == dof.get_triangulation().n_cells(level),
1690  ExcDimensionMismatch(cells.size(),
1691  dof.get_triangulation().n_cells(level)));
1692  Assert(new_order.size() == dof.n_dofs(level),
1693  ExcDimensionMismatch(new_order.size(), dof.n_dofs(level)));
1694  Assert(reverse.size() == dof.n_dofs(level),
1695  ExcDimensionMismatch(reverse.size(), dof.n_dofs(level)));
1696 
1697  unsigned int n_global_dofs = dof.n_dofs(level);
1698  unsigned int n_cell_dofs = dof.get_fe().n_dofs_per_cell();
1699 
1700  std::vector<bool> already_sorted(n_global_dofs, false);
1701  std::vector<types::global_dof_index> cell_dofs(n_cell_dofs);
1702 
1703  unsigned int global_index = 0;
1704 
1705  for (const auto &cell : cells)
1706  {
1707  Assert(cell->level() == static_cast<int>(level), ExcInternalError());
1708 
1709  cell->get_active_or_mg_dof_indices(cell_dofs);
1710  std::sort(cell_dofs.begin(), cell_dofs.end());
1711 
1712  for (unsigned int i = 0; i < n_cell_dofs; ++i)
1713  {
1714  if (!already_sorted[cell_dofs[i]])
1715  {
1716  already_sorted[cell_dofs[i]] = true;
1717  reverse[global_index++] = cell_dofs[i];
1718  }
1719  }
1720  }
1721  Assert(global_index == n_global_dofs,
1722  ExcMessage(
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?"));
1726 
1727  for (types::global_dof_index i = 0; i < new_order.size(); ++i)
1728  new_order[reverse[i]] = i;
1729  }
1730 
1731 
1732 
1733  template <typename DoFHandlerType>
1734  void
1735  downstream(DoFHandlerType & dof,
1737  const bool dof_wise_renumbering)
1738  {
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);
1743 
1744  dof.renumber_dofs(renumbering);
1745  }
1746 
1747 
1748 
1749  template <typename DoFHandlerType>
1750  void
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)
1757  {
1758  Assert(
1759  (dynamic_cast<
1760  const parallel::TriangulationBase<DoFHandlerType::dimension,
1761  DoFHandlerType::space_dimension> *>(
1762  &dof.get_triangulation()) == nullptr),
1763  ExcNotImplemented());
1764 
1765  if (dof_wise_renumbering == false)
1766  {
1767  std::vector<typename DoFHandlerType::active_cell_iterator>
1768  ordered_cells;
1769  ordered_cells.reserve(dof.get_triangulation().n_active_cells());
1770  const CompareDownstream<typename DoFHandlerType::active_cell_iterator,
1771  DoFHandlerType::space_dimension>
1772  comparator(direction);
1773 
1774  typename DoFHandlerType::active_cell_iterator p = dof.begin_active();
1775  typename DoFHandlerType::active_cell_iterator end = dof.end();
1776 
1777  while (p != end)
1778  {
1779  ordered_cells.push_back(p);
1780  ++p;
1781  }
1782  std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
1783 
1784  compute_cell_wise(new_indices, reverse, dof, ordered_cells);
1785  }
1786  else
1787  {
1788  // similar code as for
1789  // DoFTools::map_dofs_to_support_points, but
1790  // need to do this for general DoFHandlerType classes and
1791  // want to be able to sort the result
1792  // (otherwise, could use something like
1793  // DoFTools::map_support_points_to_dofs)
1794  const unsigned int n_dofs = dof.n_dofs();
1795  std::vector<
1796  std::pair<Point<DoFHandlerType::space_dimension>, unsigned int>>
1797  support_point_list(n_dofs);
1798 
1799  const hp::FECollection<DoFHandlerType::dimension> &fe_collection =
1800  dof.get_fe_collection();
1801  Assert(fe_collection[0].has_support_points(),
1802  typename FiniteElement<
1803  DoFHandlerType::dimension>::ExcFEHasNoSupportPoints());
1804  hp::QCollection<DoFHandlerType::dimension> quadrature_collection;
1805  for (unsigned int comp = 0; comp < fe_collection.size(); ++comp)
1806  {
1807  Assert(fe_collection[comp].has_support_points(),
1808  typename FiniteElement<
1809  DoFHandlerType::dimension>::ExcFEHasNoSupportPoints());
1810  quadrature_collection.push_back(
1812  fe_collection[comp].get_unit_support_points()));
1813  }
1815  hp_fe_values(fe_collection,
1816  quadrature_collection,
1818 
1819  std::vector<bool> already_touched(n_dofs, false);
1820 
1821  std::vector<types::global_dof_index> local_dof_indices;
1822  typename DoFHandlerType::active_cell_iterator begin =
1823  dof.begin_active();
1824  typename DoFHandlerType::active_cell_iterator end = dof.end();
1825  for (; begin != end; ++begin)
1826  {
1827  const unsigned int dofs_per_cell = begin->get_fe().dofs_per_cell;
1828  local_dof_indices.resize(dofs_per_cell);
1829  hp_fe_values.reinit(begin);
1830  const FEValues<DoFHandlerType::dimension> &fe_values =
1831  hp_fe_values.get_present_fe_values();
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]])
1837  {
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;
1842  }
1843  }
1844 
1846  direction);
1847  std::sort(support_point_list.begin(),
1848  support_point_list.end(),
1849  comparator);
1850  for (types::global_dof_index i = 0; i < n_dofs; ++i)
1851  new_indices[support_point_list[i].second] = i;
1852  }
1853  }
1854 
1855 
1856 
1857  template <typename DoFHandlerType>
1858  void
1859  downstream(DoFHandlerType & dof,
1860  const unsigned int level,
1862  const bool dof_wise_renumbering)
1863  {
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);
1868 
1869  dof.renumber_dofs(level, renumbering);
1870  }
1871 
1872 
1873 
1874  template <typename DoFHandlerType>
1875  void
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)
1883  {
1884  if (dof_wise_renumbering == false)
1885  {
1886  std::vector<typename DoFHandlerType::level_cell_iterator> ordered_cells;
1887  ordered_cells.reserve(dof.get_triangulation().n_cells(level));
1888  const CompareDownstream<typename DoFHandlerType::level_cell_iterator,
1889  DoFHandlerType::space_dimension>
1890  comparator(direction);
1891 
1892  typename DoFHandlerType::level_cell_iterator p = dof.begin(level);
1893  typename DoFHandlerType::level_cell_iterator end = dof.end(level);
1894 
1895  while (p != end)
1896  {
1897  ordered_cells.push_back(p);
1898  ++p;
1899  }
1900  std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
1901 
1902  compute_cell_wise(new_indices, reverse, dof, level, ordered_cells);
1903  }
1904  else
1905  {
1906  Assert(dof.get_fe().has_support_points(),
1907  typename FiniteElement<
1908  DoFHandlerType::dimension>::ExcFEHasNoSupportPoints());
1909  const unsigned int n_dofs = dof.n_dofs(level);
1910  std::vector<
1911  std::pair<Point<DoFHandlerType::space_dimension>, unsigned int>>
1912  support_point_list(n_dofs);
1913 
1915  dof.get_fe().get_unit_support_points());
1917  fe_values(dof.get_fe(), q_dummy, update_quadrature_points);
1918 
1919  std::vector<bool> already_touched(dof.n_dofs(), false);
1920 
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);
1925  for (; begin != end; ++begin)
1926  {
1927  const typename Triangulation<
1928  DoFHandlerType::dimension,
1929  DoFHandlerType::space_dimension>::cell_iterator &begin_tria =
1930  begin;
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]])
1937  {
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;
1942  }
1943  }
1944 
1946  direction);
1947  std::sort(support_point_list.begin(),
1948  support_point_list.end(),
1949  comparator);
1950  for (types::global_dof_index i = 0; i < n_dofs; ++i)
1951  new_indices[support_point_list[i].second] = i;
1952  }
1953  }
1954 
1955 
1956 
1960  namespace internal
1961  {
1962  template <int dim>
1963  struct ClockCells
1964  {
1972  bool counter;
1973 
1977  ClockCells(const Point<dim> &center, bool counter)
1978  : center(center)
1979  , counter(counter)
1980  {}
1981 
1985  template <class DHCellIterator>
1986  bool
1987  operator()(const DHCellIterator &c1, const DHCellIterator &c2) const
1988  {
1989  // dispatch to
1990  // dimension-dependent functions
1991  return compare(c1, c2, std::integral_constant<int, dim>());
1992  }
1993 
1994  private:
1998  template <class DHCellIterator, int xdim>
1999  bool
2000  compare(const DHCellIterator &c1,
2001  const DHCellIterator &c2,
2002  std::integral_constant<int, xdim>) const
2003  {
2004  const Tensor<1, dim> v1 = c1->center() - center;
2005  const Tensor<1, dim> v2 = c2->center() - center;
2006  const double s1 = std::atan2(v1[0], v1[1]);
2007  const double s2 = std::atan2(v2[0], v2[1]);
2008  return (counter ? (s1 > s2) : (s2 > s1));
2009  }
2010 
2011 
2016  template <class DHCellIterator>
2017  bool
2018  compare(const DHCellIterator &,
2019  const DHCellIterator &,
2020  std::integral_constant<int, 1>) const
2021  {
2022  Assert(dim >= 2,
2023  ExcMessage("This operation only makes sense for dim>=2."));
2024  return false;
2025  }
2026  };
2027  } // namespace internal
2028 
2029 
2030 
2031  template <typename DoFHandlerType>
2032  void
2033  clockwise_dg(DoFHandlerType & dof,
2035  const bool counter)
2036  {
2037  std::vector<types::global_dof_index> renumbering(dof.n_dofs());
2038  compute_clockwise_dg(renumbering, dof, center, counter);
2039 
2040  dof.renumber_dofs(renumbering);
2041  }
2042 
2043 
2044 
2045  template <typename DoFHandlerType>
2046  void
2047  compute_clockwise_dg(std::vector<types::global_dof_index> &new_indices,
2048  const DoFHandlerType & dof,
2050  const bool counter)
2051  {
2052  std::vector<typename DoFHandlerType::active_cell_iterator> ordered_cells;
2053  ordered_cells.reserve(dof.get_triangulation().n_active_cells());
2055  counter);
2056 
2057  typename DoFHandlerType::active_cell_iterator p = dof.begin_active();
2058  typename DoFHandlerType::active_cell_iterator end = dof.end();
2059 
2060  while (p != end)
2061  {
2062  ordered_cells.push_back(p);
2063  ++p;
2064  }
2065  std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
2066 
2067  std::vector<types::global_dof_index> reverse(new_indices.size());
2068  compute_cell_wise(new_indices, reverse, dof, ordered_cells);
2069  }
2070 
2071 
2072 
2073  template <typename DoFHandlerType>
2074  void
2075  clockwise_dg(DoFHandlerType & dof,
2076  const unsigned int level,
2078  const bool counter)
2079  {
2080  std::vector<typename DoFHandlerType::level_cell_iterator> ordered_cells;
2081  ordered_cells.reserve(dof.get_triangulation().n_active_cells());
2083  counter);
2084 
2085  typename DoFHandlerType::level_cell_iterator p = dof.begin(level);
2086  typename DoFHandlerType::level_cell_iterator end = dof.end(level);
2087 
2088  while (p != end)
2089  {
2090  ordered_cells.push_back(p);
2091  ++p;
2092  }
2093  std::sort(ordered_cells.begin(), ordered_cells.end(), comparator);
2094 
2095  cell_wise(dof, level, ordered_cells);
2096  }
2097 
2098 
2099 
2100  template <typename DoFHandlerType>
2101  void
2102  random(DoFHandlerType &dof_handler)
2103  {
2104  std::vector<types::global_dof_index> renumbering(
2105  dof_handler.n_dofs(), numbers::invalid_dof_index);
2106  compute_random(renumbering, dof_handler);
2107 
2108  dof_handler.renumber_dofs(renumbering);
2109  }
2110 
2111 
2112 
2113  template <typename DoFHandlerType>
2114  void
2115  random(DoFHandlerType &dof_handler, const unsigned int level)
2116  {
2117  Assert(dof_handler.n_dofs(level) != numbers::invalid_dof_index,
2119 
2120  std::vector<types::global_dof_index> renumbering(
2121  dof_handler.locally_owned_mg_dofs(level).n_elements(),
2123 
2124  compute_random(renumbering, dof_handler, level);
2125 
2126  dof_handler.renumber_dofs(level, renumbering);
2127  }
2128 
2129 
2130 
2131  template <typename DoFHandlerType>
2132  void
2133  compute_random(std::vector<types::global_dof_index> &new_indices,
2134  const DoFHandlerType & dof_handler)
2135  {
2136  const types::global_dof_index n_dofs = dof_handler.n_dofs();
2137  Assert(new_indices.size() == n_dofs,
2138  ExcDimensionMismatch(new_indices.size(), n_dofs));
2139 
2140  std::iota(new_indices.begin(),
2141  new_indices.end(),
2143 
2144  // shuffle the elements; the following is essentially std::shuffle (which
2145  // is new in C++11) but with a boost URNG
2146  // we could use std::mt19937 here but doing so results in compiler-dependent
2147  // output
2148  ::boost::mt19937 random_number_generator;
2149  for (unsigned int i = 1; i < n_dofs; ++i)
2150  {
2151  // get a random number between 0 and i (inclusive)
2152  const unsigned int j =
2153  ::boost::random::uniform_int_distribution<>(0, i)(
2154  random_number_generator);
2155 
2156  // if possible, swap the elements
2157  if (i != j)
2158  std::swap(new_indices[i], new_indices[j]);
2159  }
2160  }
2161 
2162 
2163 
2164  template <typename DoFHandlerType>
2165  void
2166  compute_random(std::vector<types::global_dof_index> &new_indices,
2167  const DoFHandlerType & dof_handler,
2168  const unsigned int level)
2169  {
2170  const types::global_dof_index n_dofs = dof_handler.n_dofs(level);
2171  Assert(new_indices.size() == n_dofs,
2172  ExcDimensionMismatch(new_indices.size(), n_dofs));
2173 
2174  std::iota(new_indices.begin(),
2175  new_indices.end(),
2177 
2178  // shuffle the elements; the following is essentially std::shuffle (which
2179  // is new in C++11) but with a boost URNG
2180  // we could use std::mt19937 here but doing so results in
2181  // compiler-dependent output
2182  ::boost::mt19937 random_number_generator;
2183  for (unsigned int i = 1; i < n_dofs; ++i)
2184  {
2185  // get a random number between 0 and i (inclusive)
2186  const unsigned int j =
2187  ::boost::random::uniform_int_distribution<>(0, i)(
2188  random_number_generator);
2189 
2190  // if possible, swap the elements
2191  if (i != j)
2192  std::swap(new_indices[i], new_indices[j]);
2193  }
2194  }
2195 
2196 
2197 
2198  template <typename DoFHandlerType>
2199  void
2200  subdomain_wise(DoFHandlerType &dof_handler)
2201  {
2202  Assert(
2203  (!dynamic_cast<
2204  const parallel::TriangulationBase<DoFHandlerType::dimension,
2205  DoFHandlerType::space_dimension> *>(
2206  &dof_handler.get_triangulation())),
2207  ExcMessage(
2208  "Parallel triangulations are already enumerated according to their MPI process id."));
2209 
2210  std::vector<types::global_dof_index> renumbering(
2211  dof_handler.n_dofs(), numbers::invalid_dof_index);
2212  compute_subdomain_wise(renumbering, dof_handler);
2213 
2214  dof_handler.renumber_dofs(renumbering);
2215  }
2216 
2217 
2218 
2219  template <typename DoFHandlerType>
2220  void
2221  compute_subdomain_wise(std::vector<types::global_dof_index> &new_dof_indices,
2222  const DoFHandlerType & dof_handler)
2223  {
2224  const types::global_dof_index n_dofs = dof_handler.n_dofs();
2225  Assert(new_dof_indices.size() == n_dofs,
2226  ExcDimensionMismatch(new_dof_indices.size(), n_dofs));
2227 
2228  // first get the association of each dof
2229  // with a subdomain and determine the total
2230  // number of subdomain ids used
2231  std::vector<types::subdomain_id> subdomain_association(n_dofs);
2232  DoFTools::get_subdomain_association(dof_handler, subdomain_association);
2233  const unsigned int n_subdomains =
2234  *std::max_element(subdomain_association.begin(),
2235  subdomain_association.end()) +
2236  1;
2237 
2238  // then renumber the subdomains by first
2239  // looking at those belonging to subdomain
2240  // 0, then those of subdomain 1, etc. note
2241  // that the algorithm is stable, i.e. if
2242  // two dofs i,j have i<j and belong to the
2243  // same subdomain, then they will be in
2244  // this order also after reordering
2245  std::fill(new_dof_indices.begin(),
2246  new_dof_indices.end(),
2248  types::global_dof_index next_free_index = 0;
2249  for (types::subdomain_id subdomain = 0; subdomain < n_subdomains;
2250  ++subdomain)
2251  for (types::global_dof_index i = 0; i < n_dofs; ++i)
2252  if (subdomain_association[i] == subdomain)
2253  {
2254  Assert(new_dof_indices[i] == numbers::invalid_dof_index,
2255  ExcInternalError());
2256  new_dof_indices[i] = next_free_index;
2257  ++next_free_index;
2258  }
2259 
2260  // we should have numbered all dofs
2261  Assert(next_free_index == n_dofs, ExcInternalError());
2262  Assert(std::find(new_dof_indices.begin(),
2263  new_dof_indices.end(),
2264  numbers::invalid_dof_index) == new_dof_indices.end(),
2265  ExcInternalError());
2266  }
2267 
2268 } // namespace DoFRenumbering
2269 
2270 
2271 
2272 /*-------------- Explicit Instantiations -------------------------------*/
2273 #include "dof_renumbering.inst"
2274 
2275 
identity::type
T type
Definition: template_constraints.h:270
DoFRenumbering::compute_random
void compute_random(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &dof_handler)
Definition: dof_renumbering.cc:2133
DoFRenumbering::internal::ClockCells::ClockCells
ClockCells(const Point< dim > &center, bool counter)
Definition: dof_renumbering.cc:1977
hp::FEValues
Definition: fe_values.h:281
DoFRenumbering::internal::ClockCells::counter
bool counter
Definition: dof_renumbering.cc:1972
IndexSet
Definition: index_set.h:74
DoFRenumbering::subdomain_wise
void subdomain_wise(DoFHandlerType &dof_handler)
Definition: dof_renumbering.cc:2200
DoFHandler::locally_owned_mg_dofs
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
IndexSet::subtract_set
void subtract_set(const IndexSet &other)
Definition: index_set.cc:258
dynamic_sparsity_pattern.h
numbers::invalid_dof_index
const types::global_dof_index invalid_dof_index
Definition: types.h:206
DynamicSparsityPattern::add_entries
void add_entries(const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_unique_and_sorted=false)
Definition: dynamic_sparsity_pattern.h:1051
LinearAlgebraDealII::Vector
Vector< double > Vector
Definition: generic_linear_algebra.h:43
update_quadrature_points
@ update_quadrature_points
Transformed quadrature points.
Definition: fe_update_flags.h:122
DoFHandler::active_cell_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:329
DoFRenumbering::random
void random(DoFHandlerType &dof_handler)
Definition: dof_renumbering.cc:2102
hp::FEValuesBase< dim, dim, ::FEValues< dim, dim > >::get_present_fe_values
const ::FEValues< dim, dim > & get_present_fe_values() const
Definition: fe_values.h:615
DoFRenumbering::compute_downstream
void compute_downstream(std::vector< types::global_dof_index > &new_dof_indices, std::vector< types::global_dof_index > &reverse, const DoFHandlerType &dof_handler, const Tensor< 1, DoFHandlerType::space_dimension > &direction, const bool dof_wise_renumbering)
Definition: dof_renumbering.cc:1751
DoFHandler::level_cell_iterator
typename LevelSelector::cell_iterator level_cell_iterator
Definition: dof_handler.h:393
DoFRenumbering::boost::king_ordering
void king_ordering(DoFHandlerType &dof_handler, const bool reversed_numbering=false, const bool use_constraints=false)
Definition: dof_renumbering.cc:190
DEAL_II_DOF_INDEX_MPI_TYPE
#define DEAL_II_DOF_INDEX_MPI_TYPE
Definition: types.h:86
DoFRenumbering::internal::ClockCells::center
const Point< dim > & center
Definition: dof_renumbering.cc:1968
types.h
hp::FECollection::size
unsigned int size() const
Definition: fe_collection.h:853
SparsityTools::reorder_Cuthill_McKee
void reorder_Cuthill_McKee(const DynamicSparsityPattern &sparsity, std::vector< DynamicSparsityPattern::size_type > &new_indices, const std::vector< DynamicSparsityPattern::size_type > &starting_indices=std::vector< DynamicSparsityPattern::size_type >())
Definition: sparsity_tools.cc:580
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
Triangulation
Definition: tria.h:1109
tria.h
hp::DoFHandler::active_cell_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:312
DoFRenumbering::internal::ClockCells
Definition: dof_renumbering.cc:1963
IndexSet::size
size_type size() const
Definition: index_set.h:1635
DoFRenumbering::downstream
void downstream(DoFHandlerType &dof_handler, const Tensor< 1, DoFHandlerType::space_dimension > &direction, const bool dof_wise_renumbering=false)
Definition: dof_renumbering.cc:1735
utilities.h
hp::FECollection::n_components
unsigned int n_components() const
Definition: fe_collection.h:861
tria_iterator.h
DoFRenumbering::boost::boosttypes::size_type
graph_traits< Graph >::vertices_size_type size_type
Definition: dof_renumbering.cc:83
FiniteElement::is_primitive
bool is_primitive() const
Definition: fe.h:3273
DoFRenumbering::compute_clockwise_dg
void compute_clockwise_dg(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &dof_handler, const Point< DoFHandlerType::space_dimension > &center, const bool counter)
Definition: dof_renumbering.cc:2047
DoFRenumbering::boost::compute_minimum_degree
void compute_minimum_degree(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &, const bool reversed_numbering=false, const bool use_constraints=false)
Definition: dof_renumbering.cc:266
DoFRenumbering::internal::ClockCells::operator()
bool operator()(const DHCellIterator &c1, const DHCellIterator &c2) const
Definition: dof_renumbering.cc:1987
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
Patterns::Tools::to_string
std::string to_string(const T &t)
Definition: patterns.h:2360
hp::QCollection
Definition: q_collection.h:48
MGTools::make_sparsity_pattern
void make_sparsity_pattern(const DoFHandlerType &dof_handler, SparsityPatternType &sparsity, const unsigned int level)
Definition: mg_tools.cc:565
boost
Definition: bounding_box.h:28
DoFRenumbering::compute_component_wise
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)
Definition: dof_renumbering.cc:702
sparsity_tools.h
DoFRenumbering::boost::Cuthill_McKee
void Cuthill_McKee(DoFHandlerType &dof_handler, const bool reversed_numbering=false, const bool use_constraints=false)
Definition: dof_renumbering.cc:128
hp::DoFHandler::n_dofs
types::global_dof_index n_dofs() const
DoFRenumbering::hierarchical
void hierarchical(DoFHandlerType &dof_handler)
Definition: dof_renumbering.cc:1350
second
Point< 2 > second
Definition: grid_out.cc:4353
DoFRenumbering::ExcDoFHandlerNotInitialized
static ::ExceptionBase & ExcDoFHandlerNotInitialized()
DoFRenumbering::boost::boosttypes::Graph
adjacency_list< vecS, vecS, undirectedS, property< vertex_color_t, default_color_type, property< vertex_degree_t, int > >> Graph
Definition: dof_renumbering.cc:81
DoFRenumbering::internal::ClockCells::compare
bool compare(const DHCellIterator &c1, const DHCellIterator &c2, std::integral_constant< int, xdim >) const
Definition: dof_renumbering.cc:2000
DoFHandler::renumber_dofs
void renumber_dofs(const std::vector< types::global_dof_index > &new_numbers)
Definition: dof_handler.cc:1363
hp::DoFHandler
Definition: dof_handler.h:203
DoFHandler
Definition: dof_handler.h:205
quadrature_lib.h
DoFRenumbering::internal::ClockCells::compare
bool compare(const DHCellIterator &, const DHCellIterator &, std::integral_constant< int, 1 >) const
Definition: dof_renumbering.cc:2018
Utilities::reverse_permutation
std::vector< Integer > reverse_permutation(const std::vector< Integer > &permutation)
Definition: utilities.h:1431
FEValues
Definition: fe.h:38
FiniteElement::get_nonzero_components
const ComponentMask & get_nonzero_components(const unsigned int i) const
Definition: fe.h:3253
IndexSet::n_elements
size_type n_elements() const
Definition: index_set.h:1833
DoFRenumbering::compute_cell_wise
void compute_cell_wise(std::vector< types::global_dof_index > &renumbering, std::vector< types::global_dof_index > &inverse_renumbering, const DoFHandlerType &dof_handler, const std::vector< typename DoFHandlerType::active_cell_iterator > &cell_order)
DoFRenumbering::sort_selected_dofs_back
void sort_selected_dofs_back(DoFHandlerType &dof_handler, const std::vector< bool > &selected_dofs)
Definition: dof_renumbering.cc:1455
level
unsigned int level
Definition: grid_out.cc:4355
Differentiation::SD::atan2
Expression atan2(const Expression &y, const Expression &x)
Definition: symengine_math.cc:154
DoFRenumbering::compute_block_wise
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)
Definition: dof_renumbering.cc:1058
DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition: config.h:409
DynamicSparsityPattern::row_length
size_type row_length(const size_type row) const
Definition: dynamic_sparsity_pattern.h:1072
DoFTools::make_sparsity_pattern
void make_sparsity_pattern(const DoFHandlerType &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)
Definition: dof_tools_sparsity.cc:63
DoFHandler::begin
cell_iterator begin(const unsigned int level=0) const
Definition: dof_handler.cc:922
mg_tools.h
DoFRenumbering::boost::compute_king_ordering
void compute_king_ordering(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &, const bool reversed_numbering=false, const bool use_constraints=false)
Definition: dof_renumbering.cc:210
dof_renumbering.h
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
TrilinosWrappers::internal::begin
VectorType::value_type * begin(VectorType &V)
Definition: trilinos_sparse_matrix.cc:51
DoFTools::get_subdomain_association
void get_subdomain_association(const DoFHandlerType &dof_handler, std::vector< types::subdomain_id > &subdomain)
Definition: dof_tools.cc:1550
fe_collection.h
FiniteElement::system_to_component_index
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
Definition: fe.h:3082
DoFRenumbering::boost::boosttypes::Pair
std::pair< size_type, size_type > Pair
Definition: dof_renumbering.cc:85
FiniteElementData::dofs_per_cell
const unsigned int dofs_per_cell
Definition: fe_base.h:282
hp::DoFHandler::end
cell_iterator end() const
Definition: dof_handler.cc:1338
Tensor
Definition: tensor.h:450
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
DoFRenumbering::CompareDownstream
Definition: dof_renumbering.h:353
DynamicSparsityPattern::n_rows
size_type n_rows() const
Definition: dynamic_sparsity_pattern.h:1016
FEValues::reinit
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell)
fe.h
parallel::distributed::Triangulation
Definition: tria.h:242
DynamicSparsityPattern
Definition: dynamic_sparsity_pattern.h:323
hp::FECollection::n_blocks
unsigned int n_blocks() const
Definition: fe_collection.cc:531
IndexSet::is_element
bool is_element(const size_type index) const
Definition: index_set.h:1766
hp::FEValues::reinit
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType, 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)
Definition: fe_values.cc:229
TrilinosWrappers::internal::end
VectorType::value_type * end(VectorType &V)
Definition: trilinos_sparse_matrix.cc:65
Utilities::MPI::sum
T sum(const T &t, const MPI_Comm &mpi_communicator)
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
DoFRenumbering::block_wise
void block_wise(DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_renumbering.cc:951
sparsity_pattern.h
MemorySpace::swap
void swap(MemorySpaceData< Number, MemorySpace > &, MemorySpaceData< Number, MemorySpace > &)
Definition: memory_space.h:103
AssertThrowMPI
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1707
DoFTools::make_hanging_node_constraints
void make_hanging_node_constraints(const DoFHandlerType &dof_handler, AffineConstraints< number > &constraints)
Definition: dof_tools_constraints.cc:1725
ComponentMask::first_selected_component
unsigned int first_selected_component(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
FiniteElement
Definition: fe.h:648
AffineConstraints::reinit
void reinit(const IndexSet &local_constraints=IndexSet())
DoFTools::extract_locally_active_dofs
void extract_locally_active_dofs(const DoFHandlerType &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:1138
DoFRenumbering::compute_subdomain_wise
void compute_subdomain_wise(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &dof_handler)
Definition: dof_renumbering.cc:2221
DoFRenumbering
Definition: dof_renumbering.h:345
hp::DoFHandler::renumber_dofs
void renumber_dofs(const std::vector< types::global_dof_index > &new_numbers)
Definition: dof_handler.cc:1805
DoFHandler::end
cell_iterator end() const
Definition: dof_handler.cc:951
hp::DoFHandler::begin_active
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: dof_handler.cc:1322
unsigned int
vertices
Point< 3 > vertices[4]
Definition: data_out_base.cc:174
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
DoFRenumbering::clockwise_dg
void clockwise_dg(DoFHandlerType &dof_handler, const Point< DoFHandlerType::space_dimension > &center, const bool counter=false)
Definition: dof_renumbering.cc:2033
tria.h
AffineConstraints< double >
dof_tools.h
DoFRenumbering::cell_wise
void cell_wise(DoFHandlerType &dof_handler, const std::vector< typename DoFHandlerType::active_cell_iterator > &cell_order)
Definition: dof_renumbering.cc:1570
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
hp::DoFHandler::level_cell_iterator
typename LevelSelector::cell_iterator level_cell_iterator
Definition: dof_handler.h:318
IndexSet::nth_index_in_set
size_type nth_index_in_set(const size_type local_index) const
Definition: index_set.h:1881
DoFTools::extract_locally_relevant_dofs
void extract_locally_relevant_dofs(const DoFHandlerType &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:1173
dof_handler.h
dof_accessor.h
affine_constraints.h
StandardExceptions::ExcDimensionMismatch
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Point< dim >
IndexSet::index_within_set
size_type index_within_set(const size_type global_index) const
Definition: index_set.h:1922
parallel::TriangulationBase
Definition: tria_base.h:78
DoFRenumbering::boost::minimum_degree
void minimum_degree(DoFHandlerType &dof_handler, const bool reversed_numbering=false, const bool use_constraints=false)
Definition: dof_renumbering.cc:246
DynamicSparsityPattern::column_number
size_type column_number(const size_type row, const size_type index) const
Definition: dynamic_sparsity_pattern.h:1090
template_constraints.h
FiniteElement::system_to_block_index
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
Definition: fe.h:3225
DoFHandler::n_locally_owned_dofs
types::global_dof_index n_locally_owned_dofs() const
DoFRenumbering::ComparePointwiseDownstream
Definition: dof_renumbering.h:387
internal
Definition: aligned_vector.h:369
DoFRenumbering::compute_sort_selected_dofs_back
void compute_sort_selected_dofs_back(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &dof_handler, const std::vector< bool > &selected_dofs)
Definition: dof_renumbering.cc:1490
fe_values.h
DoFRenumbering::boost::boosttypes::Vertex
graph_traits< Graph >::vertex_descriptor Vertex
Definition: dof_renumbering.cc:82
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Quadrature
Definition: quadrature.h:85
DoFHandler::begin_active
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: dof_handler.cc:935
DoFRenumbering::component_wise
void component_wise(DoFHandlerType &dof_handler, const std::vector< unsigned int > &target_component=std::vector< unsigned int >())
Definition: dof_renumbering.cc:633
TrilinosWrappers::global_index
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
Definition: trilinos_index_access.h:82
AffineConstraints::close
void close()
center
Point< 3 > center
Definition: data_out_base.cc:169
DoFRenumbering::boost::compute_Cuthill_McKee
void compute_Cuthill_McKee(std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &, const bool reversed_numbering=false, const bool use_constraints=false)
Definition: dof_renumbering.cc:148
hp::FECollection
Definition: fe_collection.h:54
internal::VectorOperations::copy
void copy(const T *begin, const T *end, U *dest)
Definition: vector_operations_internal.h:67
DoFRenumbering::boost::internal::create_graph
void create_graph(const DoFHandlerType &dof_handler, const bool use_constraints, boosttypes::Graph &graph, boosttypes::property_map< boosttypes::Graph, boosttypes::vertex_degree_t >::type &graph_degree)
Definition: dof_renumbering.cc:93
DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition: config.h:372
DoFHandler::n_dofs
types::global_dof_index n_dofs() const
hp::QCollection::push_back
void push_back(const Quadrature< dim > &new_quadrature)
Definition: q_collection.h:245
dof_handler.h