Reference documentation for deal.II version GIT 741a3088b8 2023-04-01 07:05:02+00:00
\(\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_tools.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2022 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 
16 #include <deal.II/base/mpi.h>
18 #include <deal.II/base/table.h>
20 
23 
26 #include <deal.II/dofs/dof_tools.h>
27 
28 #include <deal.II/fe/fe.h>
29 #include <deal.II/fe/fe_tools.h>
30 #include <deal.II/fe/fe_values.h>
31 
34 #include <deal.II/grid/tria.h>
36 
38 #include <deal.II/hp/fe_values.h>
41 
45 #include <deal.II/lac/vector.h>
46 
47 #include <algorithm>
48 #include <numeric>
49 
51 
52 
53 
54 namespace DoFTools
55 {
56  namespace internal
57  {
72  template <int dim, typename Number = double>
74  {
80  bool
82  const Point<dim, Number> &rhs) const
83  {
84  double downstream_size = 0;
85  double weight = 1.;
86  for (unsigned int d = 0; d < dim; ++d)
87  {
88  downstream_size += (rhs[d] - lhs[d]) * weight;
89  weight *= 1e-5;
90  }
91  if (downstream_size < 0)
92  return false;
93  else if (downstream_size > 0)
94  return true;
95  else
96  {
97  for (unsigned int d = 0; d < dim; ++d)
98  {
99  if (lhs[d] == rhs[d])
100  continue;
101  return lhs[d] < rhs[d];
102  }
103  return false;
104  }
105  }
106  };
107 
108 
109 
110  // return an array that for each dof on the reference cell
111  // lists the corresponding vector component.
112  //
113  // if an element is non-primitive then we assign to each degree of freedom
114  // the following component:
115  // - if the nonzero components that belong to a shape function are not
116  // selected in the component_mask, then the shape function is assigned
117  // to the first nonzero vector component that corresponds to this
118  // shape function
119  // - otherwise, the shape function is assigned the first component selected
120  // in the component_mask that corresponds to this shape function
121  template <int dim, int spacedim>
122  std::vector<unsigned char>
124  const ComponentMask &component_mask)
125  {
126  std::vector<unsigned char> local_component_association(
127  fe.n_dofs_per_cell(), static_cast<unsigned char>(-1));
128 
129  // compute the component each local dof belongs to.
130  // if the shape function is primitive, then this
131  // is simple and we can just associate it with
132  // what system_to_component_index gives us
133  for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
134  if (fe.is_primitive(i))
135  local_component_association[i] =
136  fe.system_to_component_index(i).first;
137  else
138  // if the shape function is not primitive, then either use the
139  // component of the first nonzero component corresponding
140  // to this shape function (if the component is not specified
141  // in the component_mask), or the first component of this block
142  // that is listed in the component_mask (if the block this
143  // component corresponds to is indeed specified in the component
144  // mask)
145  {
146  const unsigned int first_comp =
148 
149  if ((fe.get_nonzero_components(i) & component_mask)
150  .n_selected_components(fe.n_components()) == 0)
151  local_component_association[i] = first_comp;
152  else
153  // pick the component selected. we know from the previous 'if'
154  // that within the components that are nonzero for this
155  // shape function there must be at least one for which the
156  // mask is true, so we will for sure run into the break()
157  // at one point
158  for (unsigned int c = first_comp; c < fe.n_components(); ++c)
159  if (component_mask[c] == true)
160  {
161  local_component_association[i] = c;
162  break;
163  }
164  }
165 
166  Assert(std::find(local_component_association.begin(),
167  local_component_association.end(),
168  static_cast<unsigned char>(-1)) ==
169  local_component_association.end(),
170  ExcInternalError());
171 
172  return local_component_association;
173  }
174 
175 
176  // this internal function assigns to each dof the respective component
177  // of the vector system.
178  //
179  // the output array dofs_by_component lists for each dof the
180  // corresponding vector component. if the DoFHandler is based on a
181  // parallel distributed triangulation then the output array is index by
182  // dof.locally_owned_dofs().index_within_set(indices[i])
183  //
184  // if an element is non-primitive then we assign to each degree of
185  // freedom the following component:
186  // - if the nonzero components that belong to a shape function are not
187  // selected in the component_mask, then the shape function is assigned
188  // to the first nonzero vector component that corresponds to this
189  // shape function
190  // - otherwise, the shape function is assigned the first component selected
191  // in the component_mask that corresponds to this shape function
192  template <int dim, int spacedim>
193  void
195  const ComponentMask & component_mask,
196  std::vector<unsigned char> &dofs_by_component)
197  {
198  const ::hp::FECollection<dim, spacedim> &fe_collection =
199  dof.get_fe_collection();
200  Assert(fe_collection.n_components() < 256, ExcNotImplemented());
201  Assert(dofs_by_component.size() == dof.n_locally_owned_dofs(),
202  ExcDimensionMismatch(dofs_by_component.size(),
203  dof.n_locally_owned_dofs()));
204 
205  // next set up a table for the degrees of freedom on each of the
206  // cells (regardless of the fact whether it is listed in the
207  // component_select argument or not)
208  //
209  // for each element 'f' of the FECollection,
210  // local_component_association[f][d] then returns the vector
211  // component that degree of freedom 'd' belongs to
212  std::vector<std::vector<unsigned char>> local_component_association(
213  fe_collection.size());
214  for (unsigned int f = 0; f < fe_collection.size(); ++f)
215  {
216  const FiniteElement<dim, spacedim> &fe = fe_collection[f];
217  local_component_association[f] =
218  get_local_component_association(fe, component_mask);
219  }
220 
221  // then loop over all cells and do the work
222  std::vector<types::global_dof_index> indices;
223  for (const auto &c :
225  {
226  const types::fe_index fe_index = c->active_fe_index();
227  const unsigned int dofs_per_cell = c->get_fe().n_dofs_per_cell();
228  indices.resize(dofs_per_cell);
229  c->get_dof_indices(indices);
230  for (unsigned int i = 0; i < dofs_per_cell; ++i)
231  if (dof.locally_owned_dofs().is_element(indices[i]))
232  dofs_by_component[dof.locally_owned_dofs().index_within_set(
233  indices[i])] = local_component_association[fe_index][i];
234  }
235  }
236 
237 
238  // this is the function corresponding to the one above but working on
239  // blocks instead of components.
240  //
241  // the output array dofs_by_block lists for each dof the corresponding
242  // vector block. if the DoFHandler is based on a parallel distributed
243  // triangulation then the output array is index by
244  // dof.locally_owned_dofs().index_within_set(indices[i])
245  template <int dim, int spacedim>
246  inline void
248  std::vector<unsigned char> & dofs_by_block)
249  {
250  const ::hp::FECollection<dim, spacedim> &fe_collection =
251  dof.get_fe_collection();
252  Assert(fe_collection.n_components() < 256, ExcNotImplemented());
253  Assert(dofs_by_block.size() == dof.n_locally_owned_dofs(),
254  ExcDimensionMismatch(dofs_by_block.size(),
255  dof.n_locally_owned_dofs()));
256 
257  // next set up a table for the degrees of freedom on each of the
258  // cells (regardless of the fact whether it is listed in the
259  // component_select argument or not)
260  //
261  // for each element 'f' of the FECollection,
262  // local_block_association[f][d] then returns the vector block that
263  // degree of freedom 'd' belongs to
264  std::vector<std::vector<unsigned char>> local_block_association(
265  fe_collection.size());
266  for (unsigned int f = 0; f < fe_collection.size(); ++f)
267  {
268  const FiniteElement<dim, spacedim> &fe = fe_collection[f];
269  local_block_association[f].resize(fe.n_dofs_per_cell(),
270  static_cast<unsigned char>(-1));
271  for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
272  local_block_association[f][i] = fe.system_to_block_index(i).first;
273 
274  Assert(std::find(local_block_association[f].begin(),
275  local_block_association[f].end(),
276  static_cast<unsigned char>(-1)) ==
277  local_block_association[f].end(),
278  ExcInternalError());
279  }
280 
281  // then loop over all cells and do the work
282  std::vector<types::global_dof_index> indices;
283  for (const auto &cell : dof.active_cell_iterators())
284  if (cell->is_locally_owned())
285  {
286  const types::fe_index fe_index = cell->active_fe_index();
287  const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
288  indices.resize(dofs_per_cell);
289  cell->get_dof_indices(indices);
290  for (unsigned int i = 0; i < dofs_per_cell; ++i)
291  if (dof.locally_owned_dofs().is_element(indices[i]))
292  dofs_by_block[dof.locally_owned_dofs().index_within_set(
293  indices[i])] = local_block_association[fe_index][i];
294  }
295  }
296  } // namespace internal
297 
298 
299 
300  template <int dim, int spacedim, typename Number>
301  void
303  const Vector<Number> & cell_data,
304  Vector<double> & dof_data,
305  const unsigned int component)
306  {
307  const Triangulation<dim, spacedim> &tria = dof_handler.get_triangulation();
308  (void)tria;
309 
310  AssertDimension(cell_data.size(), tria.n_active_cells());
311  AssertDimension(dof_data.size(), dof_handler.n_dofs());
312  const auto &fe_collection = dof_handler.get_fe_collection();
313  AssertIndexRange(component, fe_collection.n_components());
314  for (unsigned int i = 0; i < fe_collection.size(); ++i)
315  {
316  Assert(fe_collection[i].is_primitive() == true,
318  }
319 
320  // store a flag whether we should care about different components. this
321  // is just a simplification, we could ask for this at every single
322  // place equally well
323  const bool consider_components =
324  (dof_handler.get_fe_collection().n_components() != 1);
325 
326  // zero out the components that we will touch
327  if (consider_components == false)
328  dof_data = 0;
329  else
330  {
331  std::vector<unsigned char> component_dofs(
332  dof_handler.n_locally_owned_dofs());
334  dof_handler,
335  fe_collection.component_mask(FEValuesExtractors::Scalar(component)),
336  component_dofs);
337 
338  for (unsigned int i = 0; i < dof_data.size(); ++i)
339  if (component_dofs[i] == static_cast<unsigned char>(component))
340  dof_data(i) = 0;
341  }
342 
343  // count how often we have added a value in the sum for each dof
344  std::vector<unsigned char> touch_count(dof_handler.n_dofs(), 0);
345 
346  std::vector<types::global_dof_index> dof_indices;
347  dof_indices.reserve(fe_collection.max_dofs_per_cell());
348 
349  for (const auto &cell : dof_handler.active_cell_iterators())
350  {
351  const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
352  dof_indices.resize(dofs_per_cell);
353  cell->get_dof_indices(dof_indices);
354 
355  for (unsigned int i = 0; i < dofs_per_cell; ++i)
356  // consider this dof only if it is the right component. if there
357  // is only one component, short cut the test
358  if (!consider_components ||
359  (cell->get_fe().system_to_component_index(i).first == component))
360  {
361  // sum up contribution of the present_cell to this dof
362  dof_data(dof_indices[i]) += cell_data(cell->active_cell_index());
363 
364  // note that we added another summand
365  ++touch_count[dof_indices[i]];
366  }
367  }
368 
369  // compute the mean value on all the dofs by dividing with the number
370  // of summands.
371  for (types::global_dof_index i = 0; i < dof_handler.n_dofs(); ++i)
372  {
373  // assert that each dof was used at least once. this needs not be
374  // the case if the vector has more than one component
375  Assert(consider_components || (touch_count[i] != 0),
376  ExcInternalError());
377  if (touch_count[i] != 0)
378  dof_data(i) /= touch_count[i];
379  }
380  }
381 
382 
383 
384  template <int dim, int spacedim>
385  IndexSet
387  const ComponentMask & component_mask)
388  {
389  Assert(component_mask.represents_n_components(
391  ExcMessage(
392  "The given component mask is not sized correctly to represent the "
393  "components of the given finite element."));
394 
395  // Two special cases: no component is selected, and all components are
396  // selected; both rather stupid, but easy to catch
397  if (component_mask.n_selected_components(
398  dof.get_fe_collection().n_components()) == 0)
399  return IndexSet(dof.n_dofs());
400  else if (component_mask.n_selected_components(
401  dof.get_fe_collection().n_components()) ==
403  return dof.locally_owned_dofs();
404 
405  // get the component association of each DoF and then select the ones
406  // that match the given set of components
407  std::vector<unsigned char> dofs_by_component(dof.n_locally_owned_dofs());
408  internal::get_component_association(dof, component_mask, dofs_by_component);
409 
410  // fill the selected components in a vector
411  std::vector<types::global_dof_index> selected_dofs;
412  selected_dofs.reserve(dof.n_locally_owned_dofs());
413  for (types::global_dof_index i = 0; i < dofs_by_component.size(); ++i)
414  if (component_mask[dofs_by_component[i]] == true)
415  selected_dofs.push_back(dof.locally_owned_dofs().nth_index_in_set(i));
416 
417  // fill vector of indices to return argument
418  IndexSet result(dof.n_dofs());
419  result.add_indices(selected_dofs.begin(), selected_dofs.end());
420  return result;
421  }
422 
423 
424 
425  template <int dim, int spacedim>
426  IndexSet
428  const BlockMask & block_mask)
429  {
430  // simply forward to the function that works based on a component mask
431  return extract_dofs<dim, spacedim>(
432  dof, dof.get_fe_collection().component_mask(block_mask));
433  }
434 
435 
436 
437  template <int dim, int spacedim>
438  std::vector<IndexSet>
440  const ComponentMask &component_mask)
441  {
442  const auto n_comps = dof.get_fe_collection().n_components();
443  Assert(component_mask.represents_n_components(n_comps),
444  ExcMessage(
445  "The given component mask is not sized correctly to represent the "
446  "components of the given finite element."));
447 
448  const auto &locally_owned_dofs = dof.locally_owned_dofs();
449 
450  // get the component association of each DoF and then select the ones
451  // that match the given set of components
452  std::vector<unsigned char> dofs_by_component(dof.n_locally_owned_dofs());
453  internal::get_component_association(dof, component_mask, dofs_by_component);
454 
455  std::vector<IndexSet> index_per_comp(n_comps, IndexSet(dof.n_dofs()));
456 
457  for (types::global_dof_index i = 0; i < dof.n_locally_owned_dofs(); ++i)
458  {
459  const auto &comp_i = dofs_by_component[i];
460  if (component_mask[comp_i])
461  index_per_comp[comp_i].add_index(
462  locally_owned_dofs.nth_index_in_set(i));
463  }
464  for (auto &c : index_per_comp)
465  c.compress();
466  return index_per_comp;
467  }
468 
469 
470 
471  template <int dim, int spacedim>
472  void
473  extract_level_dofs(const unsigned int level,
474  const DoFHandler<dim, spacedim> &dof,
475  const ComponentMask & component_mask,
476  std::vector<bool> & selected_dofs)
477  {
478  const FiniteElement<dim, spacedim> &fe = dof.get_fe();
479 
480  Assert(component_mask.represents_n_components(
482  ExcMessage(
483  "The given component mask is not sized correctly to represent the "
484  "components of the given finite element."));
485  Assert(selected_dofs.size() == dof.n_dofs(level),
486  ExcDimensionMismatch(selected_dofs.size(), dof.n_dofs(level)));
487 
488  // two special cases: no component is selected, and all components are
489  // selected, both rather stupid, but easy to catch
490  if (component_mask.n_selected_components(
491  dof.get_fe_collection().n_components()) == 0)
492  {
493  std::fill_n(selected_dofs.begin(), dof.n_dofs(level), false);
494  return;
495  }
496  else if (component_mask.n_selected_components(
497  dof.get_fe_collection().n_components()) ==
499  {
500  std::fill_n(selected_dofs.begin(), dof.n_dofs(level), true);
501  return;
502  }
503 
504  // preset all values by false
505  std::fill_n(selected_dofs.begin(), dof.n_dofs(level), false);
506 
507  // next set up a table for the degrees of freedom on each of the cells
508  // whether it is something interesting or not
509  std::vector<unsigned char> local_component_asssociation =
510  internal::get_local_component_association(fe, component_mask);
511  std::vector<bool> local_selected_dofs(fe.n_dofs_per_cell());
512  for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
513  local_selected_dofs[i] = component_mask[local_component_asssociation[i]];
514 
515  // then loop over all cells and do work
516  std::vector<types::global_dof_index> indices(fe.n_dofs_per_cell());
517  for (const auto &cell : dof.cell_iterators_on_level(level))
518  {
519  cell->get_mg_dof_indices(indices);
520  for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
521  selected_dofs[indices[i]] = local_selected_dofs[i];
522  }
523  }
524 
525 
526 
527  template <int dim, int spacedim>
528  void
529  extract_level_dofs(const unsigned int level,
530  const DoFHandler<dim, spacedim> &dof,
531  const BlockMask & block_mask,
532  std::vector<bool> & selected_dofs)
533  {
534  // simply defer to the other extract_level_dofs() function
536  dof,
537  dof.get_fe().component_mask(block_mask),
538  selected_dofs);
539  }
540 
541 
542 
543  template <int dim, int spacedim>
544  void
546  const ComponentMask & component_mask,
547  std::vector<bool> & selected_dofs,
548  const std::set<types::boundary_id> &boundary_ids)
549  {
550  Assert((dynamic_cast<
552  &dof_handler.get_triangulation()) == nullptr),
553  ExcMessage(
554  "This function can not be used with distributed triangulations. "
555  "See the documentation for more information."));
556 
557  IndexSet indices;
558  extract_boundary_dofs(dof_handler, component_mask, indices, boundary_ids);
559 
560  // clear and reset array by default values
561  selected_dofs.clear();
562  selected_dofs.resize(dof_handler.n_dofs(), false);
563 
564  // then convert the values computed above to the binary vector
565  indices.fill_binary_vector(selected_dofs);
566  }
567 
568 
569 
570  template <int dim, int spacedim>
571  void
573  const ComponentMask & component_mask,
574  IndexSet & selected_dofs,
575  const std::set<types::boundary_id> &boundary_ids)
576  {
577  // Simply forward to the other function
578  selected_dofs =
579  extract_boundary_dofs(dof_handler, component_mask, boundary_ids);
580  }
581 
582 
583 
584  template <int dim, int spacedim>
585  IndexSet
587  const ComponentMask & component_mask,
588  const std::set<types::boundary_id> &boundary_ids)
589  {
590  Assert(component_mask.represents_n_components(
591  dof_handler.get_fe_collection().n_components()),
592  ExcMessage("Component mask has invalid size."));
593  Assert(boundary_ids.find(numbers::internal_face_boundary_id) ==
594  boundary_ids.end(),
596 
597  IndexSet selected_dofs(dof_handler.n_dofs());
598 
599  // let's see whether we have to check for certain boundary indicators
600  // or whether we can accept all
601  const bool check_boundary_id = (boundary_ids.size() != 0);
602 
603  // also see whether we have to check whether a certain vector component
604  // is selected, or all
605  const bool check_vector_component =
606  ((component_mask.represents_the_all_selected_mask() == false) ||
607  (component_mask.n_selected_components(
608  dof_handler.get_fe_collection().n_components()) !=
609  dof_handler.get_fe_collection().n_components()));
610 
611  std::vector<types::global_dof_index> face_dof_indices;
612  face_dof_indices.reserve(
613  dof_handler.get_fe_collection().max_dofs_per_face());
614 
615  // now loop over all cells and check whether their faces are at the
616  // boundary. note that we need not take special care of single lines
617  // being at the boundary (using @p{cell->has_boundary_lines}), since we
618  // do not support boundaries of dimension dim-2, and so every isolated
619  // boundary line is also part of a boundary face which we will be
620  // visiting sooner or later
621  for (const auto &cell : dof_handler.active_cell_iterators())
622  // only work on cells that are either locally owned or at least ghost
623  // cells
624  if (cell->is_artificial() == false)
625  for (const unsigned int face : cell->face_indices())
626  if (cell->at_boundary(face))
627  if (!check_boundary_id ||
628  (boundary_ids.find(cell->face(face)->boundary_id()) !=
629  boundary_ids.end()))
630  {
631  const FiniteElement<dim, spacedim> &fe = cell->get_fe();
632 
633  const auto reference_cell = cell->reference_cell();
634 
635  const unsigned int n_vertices_per_cell =
636  reference_cell.n_vertices();
637  const unsigned int n_lines_per_cell = reference_cell.n_lines();
638  const unsigned int n_vertices_per_face =
639  reference_cell.face_reference_cell(face).n_vertices();
640  const unsigned int n_lines_per_face =
641  reference_cell.face_reference_cell(face).n_lines();
642 
643  const unsigned int dofs_per_face = fe.n_dofs_per_face(face);
644  face_dof_indices.resize(dofs_per_face);
645  cell->face(face)->get_dof_indices(face_dof_indices,
646  cell->active_fe_index());
647 
648  for (unsigned int i = 0; i < fe.n_dofs_per_face(face); ++i)
649  if (!check_vector_component)
650  selected_dofs.add_index(face_dof_indices[i]);
651  else
652  // check for component is required. somewhat tricky as
653  // usual for the case that the shape function is
654  // non-primitive, but use usual convention (see docs)
655  {
656  // first get at the cell-global number of a face dof,
657  // to ask the FE certain questions
658  const unsigned int cell_index =
659  (dim == 1 ?
660  i :
661  (dim == 2 ?
662  (i < 2 * fe.n_dofs_per_vertex() ?
663  i :
664  i + 2 * fe.n_dofs_per_vertex()) :
665  (dim == 3 ? (i < n_vertices_per_face *
666  fe.n_dofs_per_vertex() ?
667  i :
668  (i < n_vertices_per_face *
669  fe.n_dofs_per_vertex() +
670  n_lines_per_face *
671  fe.n_dofs_per_line() ?
672  (i - n_vertices_per_face *
673  fe.n_dofs_per_vertex()) +
674  n_vertices_per_cell *
675  fe.n_dofs_per_vertex() :
676  (i -
677  n_vertices_per_face *
678  fe.n_dofs_per_vertex() -
679  n_lines_per_face *
680  fe.n_dofs_per_line()) +
681  n_vertices_per_cell *
682  fe.n_dofs_per_vertex() +
683  n_lines_per_cell *
684  fe.n_dofs_per_line())) :
686  if (fe.is_primitive(cell_index))
687  {
688  if (component_mask[fe.face_system_to_component_index(
689  i, face)
690  .first] == true)
691  selected_dofs.add_index(face_dof_indices[i]);
692  }
693  else // not primitive
694  {
695  const unsigned int first_nonzero_comp =
698  Assert(first_nonzero_comp < fe.n_components(),
699  ExcInternalError());
700 
701  if (component_mask[first_nonzero_comp] == true)
702  selected_dofs.add_index(face_dof_indices[i]);
703  }
704  }
705  }
706 
707  return selected_dofs;
708  }
709 
710 
711 
712  template <int dim, int spacedim>
713  void
715  const DoFHandler<dim, spacedim> & dof_handler,
716  const ComponentMask & component_mask,
717  std::vector<bool> & selected_dofs,
718  const std::set<types::boundary_id> &boundary_ids)
719  {
720  Assert(component_mask.represents_n_components(
721  dof_handler.get_fe_collection().n_components()),
722  ExcMessage("This component mask has the wrong size."));
723  Assert(boundary_ids.find(numbers::internal_face_boundary_id) ==
724  boundary_ids.end(),
726 
727  // let's see whether we have to check for certain boundary indicators
728  // or whether we can accept all
729  const bool check_boundary_id = (boundary_ids.size() != 0);
730 
731  // also see whether we have to check whether a certain vector component
732  // is selected, or all
733  const bool check_vector_component =
734  (component_mask.represents_the_all_selected_mask() == false);
735 
736  // clear and reset array by default values
737  selected_dofs.clear();
738  selected_dofs.resize(dof_handler.n_dofs(), false);
739  std::vector<types::global_dof_index> cell_dof_indices;
740  cell_dof_indices.reserve(
741  dof_handler.get_fe_collection().max_dofs_per_cell());
742 
743  // now loop over all cells and check whether their faces are at the
744  // boundary. note that we need not take special care of single lines
745  // being at the boundary (using @p{cell->has_boundary_lines}), since we
746  // do not support boundaries of dimension dim-2, and so every isolated
747  // boundary line is also part of a boundary face which we will be
748  // visiting sooner or later
749  for (const auto &cell : dof_handler.active_cell_iterators())
750  for (const unsigned int face : cell->face_indices())
751  if (cell->at_boundary(face))
752  if (!check_boundary_id ||
753  (boundary_ids.find(cell->face(face)->boundary_id()) !=
754  boundary_ids.end()))
755  {
756  const FiniteElement<dim, spacedim> &fe = cell->get_fe();
757 
758  const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
759  cell_dof_indices.resize(dofs_per_cell);
760  cell->get_dof_indices(cell_dof_indices);
761 
762  for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
763  if (fe.has_support_on_face(i, face))
764  {
765  if (!check_vector_component)
766  selected_dofs[cell_dof_indices[i]] = true;
767  else
768  // check for component is required. somewhat tricky
769  // as usual for the case that the shape function is
770  // non-primitive, but use usual convention (see docs)
771  {
772  if (fe.is_primitive(i))
773  selected_dofs[cell_dof_indices[i]] =
774  (component_mask[fe.system_to_component_index(i)
775  .first] == true);
776  else // not primitive
777  {
778  const unsigned int first_nonzero_comp =
781  Assert(first_nonzero_comp < fe.n_components(),
782  ExcInternalError());
783 
784  selected_dofs[cell_dof_indices[i]] =
785  (component_mask[first_nonzero_comp] == true);
786  }
787  }
788  }
789  }
790  }
791 
792 
793 
794  template <int dim, int spacedim, typename number>
795  IndexSet
797  const DoFHandler<dim, spacedim> &dof_handler,
798  const std::function<
800  & predicate,
801  const AffineConstraints<number> &cm)
802  {
803  const std::function<bool(
805  predicate_local =
806  [=](
808  -> bool { return cell->is_locally_owned() && predicate(cell); };
809 
810  std::vector<types::global_dof_index> local_dof_indices;
811  local_dof_indices.reserve(
812  dof_handler.get_fe_collection().max_dofs_per_cell());
813 
814  // Get all the dofs that live on the subdomain:
815  std::set<types::global_dof_index> predicate_dofs;
816 
817  for (const auto &cell : dof_handler.active_cell_iterators())
818  if (!cell->is_artificial() && predicate(cell))
819  {
820  local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
821  cell->get_dof_indices(local_dof_indices);
822  predicate_dofs.insert(local_dof_indices.begin(),
823  local_dof_indices.end());
824  }
825 
826  // Get halo layer and accumulate its DoFs
827  std::set<types::global_dof_index> dofs_with_support_on_halo_cells;
828 
829  const std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
830  halo_cells =
831  GridTools::compute_active_cell_halo_layer(dof_handler, predicate_local);
832  for (typename std::vector<typename DoFHandler<dim, spacedim>::
833  active_cell_iterator>::const_iterator it =
834  halo_cells.begin();
835  it != halo_cells.end();
836  ++it)
837  {
838  // skip halo cells that satisfy the given predicate and thereby
839  // shall be a part of the index set on another MPI core.
840  // Those could only be ghost cells with p::d::Tria.
841  if (predicate(*it))
842  {
843  Assert((*it)->is_ghost(), ExcInternalError());
844  continue;
845  }
846 
847  const unsigned int dofs_per_cell = (*it)->get_fe().n_dofs_per_cell();
848  local_dof_indices.resize(dofs_per_cell);
849  (*it)->get_dof_indices(local_dof_indices);
850  dofs_with_support_on_halo_cells.insert(local_dof_indices.begin(),
851  local_dof_indices.end());
852  }
853 
854  // A complication coming from the fact that DoFs living on locally
855  // owned cells which satisfy predicate may participate in constraints for
856  // DoFs outside of this region.
857  if (cm.n_constraints() > 0)
858  {
859  // Remove DoFs that are part of constraints for DoFs outside
860  // of predicate. Since we will subtract halo_dofs from predicate_dofs,
861  // achieve this by extending halo_dofs with DoFs to which
862  // halo_dofs are constrained.
863  std::set<types::global_dof_index> extra_halo;
864  for (const auto dof : dofs_with_support_on_halo_cells)
865  // if halo DoF is constrained, add all DoFs to which it's constrained
866  // because after resolving constraints, the support of the DoFs that
867  // constrain the current DoF will extend to the halo cells.
868  if (const auto *line_ptr = cm.get_constraint_entries(dof))
869  {
870  const unsigned int line_size = line_ptr->size();
871  for (unsigned int j = 0; j < line_size; ++j)
872  extra_halo.insert((*line_ptr)[j].first);
873  }
874 
875  dofs_with_support_on_halo_cells.insert(extra_halo.begin(),
876  extra_halo.end());
877  }
878 
879  // Rework our std::set's into IndexSet and subtract halo DoFs from the
880  // predicate_dofs:
881  IndexSet support_set(dof_handler.n_dofs());
882  support_set.add_indices(predicate_dofs.begin(), predicate_dofs.end());
883  support_set.compress();
884 
885  IndexSet halo_set(dof_handler.n_dofs());
886  halo_set.add_indices(dofs_with_support_on_halo_cells.begin(),
887  dofs_with_support_on_halo_cells.end());
888  halo_set.compress();
889 
890  support_set.subtract_set(halo_set);
891 
892  // we intentionally do not want to limit the output to locally owned DoFs.
893  return support_set;
894  }
895 
896 
897 
898  namespace internal
899  {
900  namespace
901  {
902  template <int spacedim>
903  IndexSet
905  {
906  // there are no hanging nodes in 1d
907  return IndexSet(dof_handler.n_dofs());
908  }
909 
910 
911  template <int spacedim>
912  IndexSet
914  {
915  const unsigned int dim = 2;
916 
917  IndexSet selected_dofs(dof_handler.n_dofs());
918 
919  const FiniteElement<dim, spacedim> &fe = dof_handler.get_fe();
920 
921  // this function is similar to the make_sparsity_pattern function,
922  // see there for more information
923  for (const auto &cell : dof_handler.active_cell_iterators())
924  if (!cell->is_artificial())
925  {
926  for (const unsigned int face : cell->face_indices())
927  if (cell->face(face)->has_children())
928  {
930  line = cell->face(face);
931 
932  for (unsigned int dof = 0; dof != fe.n_dofs_per_vertex();
933  ++dof)
934  selected_dofs.add_index(
935  line->child(0)->vertex_dof_index(1, dof));
936 
937  for (unsigned int child = 0; child < 2; ++child)
938  {
939  if (cell->neighbor_child_on_subface(face, child)
940  ->is_artificial())
941  continue;
942  for (unsigned int dof = 0; dof != fe.n_dofs_per_line();
943  ++dof)
944  selected_dofs.add_index(
945  line->child(child)->dof_index(dof));
946  }
947  }
948  }
949 
950  selected_dofs.compress();
951  return selected_dofs;
952  }
953 
954 
955  template <int spacedim>
956  IndexSet
958  {
959  const unsigned int dim = 3;
960 
961  IndexSet selected_dofs(dof_handler.n_dofs());
962  IndexSet unconstrained_dofs(dof_handler.n_dofs());
963 
964  const FiniteElement<dim, spacedim> &fe = dof_handler.get_fe();
965 
966  for (const auto &cell : dof_handler.active_cell_iterators())
967  if (!cell->is_artificial())
968  for (auto f : cell->face_indices())
969  {
970  const typename DoFHandler<dim, spacedim>::face_iterator face =
971  cell->face(f);
972  if (cell->face(f)->has_children())
973  {
974  for (unsigned int child = 0; child < 4; ++child)
975  if (!cell->neighbor_child_on_subface(f, child)
976  ->is_artificial())
977  {
978  // simply take all DoFs that live on this subface
979  std::vector<types::global_dof_index> ldi(
980  fe.n_dofs_per_face(f, child));
981  face->child(child)->get_dof_indices(ldi);
982  selected_dofs.add_indices(ldi.begin(), ldi.end());
983  }
984 
985  // and subtract (in the end) all the indices which a shared
986  // between this face and its subfaces
987  for (unsigned int vertex = 0; vertex < 4; ++vertex)
988  for (unsigned int dof = 0; dof != fe.n_dofs_per_vertex();
989  ++dof)
990  unconstrained_dofs.add_index(
991  face->vertex_dof_index(vertex, dof));
992  }
993  }
994  selected_dofs.subtract_set(unconstrained_dofs);
995  return selected_dofs;
996  }
997  } // namespace
998  } // namespace internal
999 
1000 
1001 
1002  template <int dim, int spacedim>
1003  IndexSet
1005  {
1006  return internal::extract_hanging_node_dofs(dof_handler);
1007  }
1008 
1009 
1010 
1011  template <int dim, int spacedim>
1012  void
1015  std::vector<bool> & selected_dofs)
1016  {
1017  Assert(selected_dofs.size() == dof_handler.n_dofs(),
1018  ExcDimensionMismatch(selected_dofs.size(), dof_handler.n_dofs()));
1019 
1020  // preset all values by false
1021  std::fill_n(selected_dofs.begin(), dof_handler.n_dofs(), false);
1022 
1023  std::vector<types::global_dof_index> local_dof_indices;
1024  local_dof_indices.reserve(
1025  dof_handler.get_fe_collection().max_dofs_per_cell());
1026 
1027  // this function is similar to the make_sparsity_pattern function, see
1028  // there for more information
1029  for (const auto &cell : dof_handler.active_cell_iterators())
1030  if (cell->subdomain_id() == subdomain_id)
1031  {
1032  const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
1033  local_dof_indices.resize(dofs_per_cell);
1034  cell->get_dof_indices(local_dof_indices);
1035  for (unsigned int i = 0; i < dofs_per_cell; ++i)
1036  selected_dofs[local_dof_indices[i]] = true;
1037  }
1038  }
1039 
1040 
1041 
1042  template <int dim, int spacedim>
1043  IndexSet
1045  {
1046  // collect all the locally owned dofs
1047  IndexSet dof_set = dof_handler.locally_owned_dofs();
1048 
1049  // add the DoF on the adjacent ghost cells to the IndexSet, cache them
1050  // in a set. need to check each dof manually because we can't be sure
1051  // that the dof range of locally_owned_dofs is really contiguous.
1052  std::vector<types::global_dof_index> dof_indices;
1053  std::set<types::global_dof_index> global_dof_indices;
1054 
1055  for (const auto &cell : dof_handler.active_cell_iterators() |
1057  {
1058  dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1059  cell->get_dof_indices(dof_indices);
1060 
1061  for (const types::global_dof_index dof_index : dof_indices)
1062  if (!dof_set.is_element(dof_index))
1063  global_dof_indices.insert(dof_index);
1064  }
1065 
1066  dof_set.add_indices(global_dof_indices.begin(), global_dof_indices.end());
1067 
1068  dof_set.compress();
1069 
1070  return dof_set;
1071  }
1072 
1073 
1074 
1075  template <int dim, int spacedim>
1076  void
1078  IndexSet & dof_set)
1079  {
1080  dof_set = extract_locally_active_dofs(dof_handler);
1081  }
1082 
1083 
1084 
1085  template <int dim, int spacedim>
1086  IndexSet
1088  const DoFHandler<dim, spacedim> &dof_handler,
1089  const unsigned int level)
1090  {
1091  // collect all the locally owned dofs
1092  IndexSet dof_set = dof_handler.locally_owned_mg_dofs(level);
1093 
1094  // add the DoF on the adjacent ghost cells to the IndexSet, cache them
1095  // in a set. need to check each dof manually because we can't be sure
1096  // that the dof range of locally_owned_dofs is really contiguous.
1097  std::vector<types::global_dof_index> dof_indices;
1098  std::set<types::global_dof_index> global_dof_indices;
1099 
1100  const auto filtered_iterators_range =
1103  for (const auto &cell : filtered_iterators_range)
1104  {
1105  dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1106  cell->get_mg_dof_indices(dof_indices);
1107 
1108  for (const types::global_dof_index dof_index : dof_indices)
1109  if (!dof_set.is_element(dof_index))
1110  global_dof_indices.insert(dof_index);
1111  }
1112 
1113  dof_set.add_indices(global_dof_indices.begin(), global_dof_indices.end());
1114 
1115  dof_set.compress();
1116 
1117  return dof_set;
1118  }
1119 
1120 
1121 
1122  template <int dim, int spacedim>
1123  void
1125  const DoFHandler<dim, spacedim> &dof_handler,
1126  IndexSet & dof_set,
1127  const unsigned int level)
1128  {
1129  dof_set = extract_locally_active_level_dofs(dof_handler, level);
1130  }
1131 
1132 
1133 
1134  template <int dim, int spacedim>
1135  IndexSet
1137  {
1138  // collect all the locally owned dofs
1139  IndexSet dof_set = dof_handler.locally_owned_dofs();
1140 
1141  // now add the DoF on the adjacent ghost cells to the IndexSet
1142 
1143  // Note: For certain meshes (in particular in 3d and with many
1144  // processors), it is really necessary to cache intermediate data. After
1145  // trying several objects such as std::set, a vector that is always kept
1146  // sorted, and a vector that is initially unsorted and sorted once at the
1147  // end, the latter has been identified to provide the best performance.
1148  // Martin Kronbichler
1149  std::vector<types::global_dof_index> dof_indices;
1150  std::vector<types::global_dof_index> dofs_on_ghosts;
1151 
1152  for (const auto &cell : dof_handler.active_cell_iterators())
1153  if (cell->is_ghost())
1154  {
1155  dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1156  cell->get_dof_indices(dof_indices);
1157  for (const auto dof_index : dof_indices)
1158  if (!dof_set.is_element(dof_index))
1159  dofs_on_ghosts.push_back(dof_index);
1160  }
1161 
1162  // sort, compress out duplicates, fill into index set
1163  std::sort(dofs_on_ghosts.begin(), dofs_on_ghosts.end());
1164  dof_set.add_indices(dofs_on_ghosts.begin(),
1165  std::unique(dofs_on_ghosts.begin(),
1166  dofs_on_ghosts.end()));
1167  dof_set.compress();
1168 
1169  return dof_set;
1170  }
1171 
1172 
1173 
1174  template <int dim, int spacedim>
1175  void
1177  IndexSet & dof_set)
1178  {
1179  dof_set = extract_locally_relevant_dofs(dof_handler);
1180  }
1181 
1182 
1183 
1184  template <int dim, int spacedim>
1185  IndexSet
1187  const DoFHandler<dim, spacedim> &dof_handler,
1188  const unsigned int level)
1189  {
1190  // collect all the locally owned dofs
1191  IndexSet dof_set = dof_handler.locally_owned_mg_dofs(level);
1192 
1193  // add the DoF on the adjacent ghost cells to the IndexSet
1194 
1195  // Note: For certain meshes (in particular in 3d and with many
1196  // processors), it is really necessary to cache intermediate data. After
1197  // trying several objects such as std::set, a vector that is always kept
1198  // sorted, and a vector that is initially unsorted and sorted once at the
1199  // end, the latter has been identified to provide the best performance.
1200  // Martin Kronbichler
1201  std::vector<types::global_dof_index> dof_indices;
1202  std::vector<types::global_dof_index> dofs_on_ghosts;
1203 
1204  for (const auto &cell : dof_handler.cell_iterators_on_level(level))
1205  {
1206  const types::subdomain_id id = cell->level_subdomain_id();
1207 
1208  // skip artificial and own cells (only look at ghost cells)
1209  if (id == dof_handler.get_triangulation().locally_owned_subdomain() ||
1211  continue;
1212 
1213  dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1214  cell->get_mg_dof_indices(dof_indices);
1215  for (const auto dof_index : dof_indices)
1216  if (!dof_set.is_element(dof_index))
1217  dofs_on_ghosts.push_back(dof_index);
1218  }
1219 
1220  // sort, compress out duplicates, fill into index set
1221  std::sort(dofs_on_ghosts.begin(), dofs_on_ghosts.end());
1222  dof_set.add_indices(dofs_on_ghosts.begin(),
1223  std::unique(dofs_on_ghosts.begin(),
1224  dofs_on_ghosts.end()));
1225 
1226  dof_set.compress();
1227 
1228  return dof_set;
1229  }
1230 
1231 
1232 
1233  template <int dim, int spacedim>
1234  void
1236  const DoFHandler<dim, spacedim> &dof_handler,
1237  const unsigned int level,
1238  IndexSet & dof_set)
1239  {
1240  dof_set = extract_locally_relevant_level_dofs(dof_handler, level);
1241  }
1242 
1243 
1244 
1245  template <int dim, int spacedim>
1246  void
1248  const ComponentMask & component_mask,
1249  std::vector<std::vector<bool>> & constant_modes)
1250  {
1251  // If there are no locally owned DoFs, return with an empty
1252  // constant_modes object:
1253  if (dof_handler.n_locally_owned_dofs() == 0)
1254  {
1255  constant_modes = std::vector<std::vector<bool>>(0);
1256  return;
1257  }
1258 
1259  const unsigned int n_components = dof_handler.get_fe(0).n_components();
1260  Assert(component_mask.represents_n_components(n_components),
1261  ExcDimensionMismatch(n_components, component_mask.size()));
1262 
1263  std::vector<unsigned char> dofs_by_component(
1264  dof_handler.n_locally_owned_dofs());
1266  component_mask,
1267  dofs_by_component);
1268  unsigned int n_selected_dofs = 0;
1269  for (unsigned int i = 0; i < n_components; ++i)
1270  if (component_mask[i] == true)
1271  n_selected_dofs +=
1272  std::count(dofs_by_component.begin(), dofs_by_component.end(), i);
1273 
1274  // Find local numbering within the selected components
1275  const IndexSet &locally_owned_dofs = dof_handler.locally_owned_dofs();
1276  std::vector<unsigned int> component_numbering(
1277  locally_owned_dofs.n_elements(), numbers::invalid_unsigned_int);
1278  for (unsigned int i = 0, count = 0; i < locally_owned_dofs.n_elements();
1279  ++i)
1280  if (component_mask[dofs_by_component[i]])
1281  component_numbering[i] = count++;
1282 
1283  // get the element constant modes and find a translation table between
1284  // index in the constant modes and the components.
1285  //
1286  // TODO: We might be able to extend this also for elements which do not
1287  // have the same constant modes, but that is messy...
1288  const ::hp::FECollection<dim, spacedim> &fe_collection =
1289  dof_handler.get_fe_collection();
1290  std::vector<Table<2, bool>> element_constant_modes;
1291  std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
1292  constant_mode_to_component_translation(n_components);
1293  {
1294  unsigned int n_constant_modes = 0;
1295  int first_non_empty_constant_mode = -1;
1296  for (unsigned int f = 0; f < fe_collection.size(); ++f)
1297  {
1298  std::pair<Table<2, bool>, std::vector<unsigned int>> data =
1299  fe_collection[f].get_constant_modes();
1300 
1301  // Store the index of the current element if it is the first that has
1302  // non-empty constant modes.
1303  if (first_non_empty_constant_mode < 0 && data.first.n_rows() > 0)
1304  {
1305  first_non_empty_constant_mode = f;
1306  // This is the first non-empty constant mode, so we figure out the
1307  // translation between index in the constant modes and the
1308  // components
1309  for (unsigned int i = 0; i < data.second.size(); ++i)
1310  if (component_mask[data.second[i]])
1311  constant_mode_to_component_translation[data.second[i]]
1312  .emplace_back(n_constant_modes++, i);
1313  }
1314 
1315  // Add the constant modes of this element to the list and assert that
1316  // there are as many constant modes as for the other elements (or zero
1317  // constant modes).
1318  element_constant_modes.push_back(data.first);
1319  Assert(
1320  element_constant_modes.back().n_rows() == 0 ||
1321  element_constant_modes.back().n_rows() ==
1322  element_constant_modes[first_non_empty_constant_mode].n_rows(),
1323  ExcInternalError());
1324  }
1325  AssertIndexRange(first_non_empty_constant_mode, fe_collection.size());
1326 
1327  // Now we know the number of constant modes and resize the return vector
1328  // accordingly
1329  constant_modes.clear();
1330  constant_modes.resize(n_constant_modes,
1331  std::vector<bool>(n_selected_dofs, false));
1332  }
1333 
1334  // Loop over all owned cells and ask the element for the constant modes
1335  std::vector<types::global_dof_index> dof_indices;
1336  for (const auto &cell : dof_handler.active_cell_iterators() |
1338  {
1339  dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1340  cell->get_dof_indices(dof_indices);
1341 
1342  for (unsigned int i = 0; i < dof_indices.size(); ++i)
1343  if (locally_owned_dofs.is_element(dof_indices[i]))
1344  {
1345  const unsigned int loc_index =
1346  locally_owned_dofs.index_within_set(dof_indices[i]);
1347  const unsigned int comp = dofs_by_component[loc_index];
1348  if (component_mask[comp])
1349  for (auto &indices :
1350  constant_mode_to_component_translation[comp])
1351  constant_modes[indices
1352  .first][component_numbering[loc_index]] =
1353  element_constant_modes[cell->active_fe_index()](
1354  indices.second, i);
1355  }
1356  }
1357  }
1358 
1359 
1360 
1361  template <int dim, int spacedim>
1362  void
1364  std::vector<unsigned int> & active_fe_indices)
1365  {
1366  AssertDimension(active_fe_indices.size(),
1367  dof_handler.get_triangulation().n_active_cells());
1368 
1369  std::vector<types::fe_index> indices = dof_handler.get_active_fe_indices();
1370 
1371  active_fe_indices.assign(indices.begin(), indices.end());
1372  }
1373 
1374  template <int dim, int spacedim>
1375  std::vector<IndexSet>
1377  {
1378  Assert(dof_handler.n_dofs() > 0,
1379  ExcMessage("The given DoFHandler has no DoFs."));
1380 
1381  // If the Triangulation is distributed, the only thing we can usefully
1382  // ask is for its locally owned subdomain
1383  Assert((dynamic_cast<
1385  &dof_handler.get_triangulation()) == nullptr),
1386  ExcMessage(
1387  "For parallel::distributed::Triangulation objects and "
1388  "associated DoF handler objects, asking for any information "
1389  "related to a subdomain other than the locally owned one does "
1390  "not make sense."));
1391 
1392  // The following is a random process (flip of a coin), thus should be called
1393  // once only.
1394  std::vector<::types::subdomain_id> subdomain_association(
1395  dof_handler.n_dofs());
1397  subdomain_association);
1398 
1399  // Figure out how many subdomain ids there are.
1400  //
1401  // if this is a parallel triangulation, then we can just ask the
1402  // triangulation for this. if this is a sequential triangulation, we loop
1403  // over all cells and take the largest subdomain_id value we find; the
1404  // number of subdomains is then the largest found value plus one. (we here
1405  // assume that all subdomain ids up to the largest are actually used; this
1406  // may not be true for a sequential triangulation where these values have
1407  // been set by hand and not in accordance with some MPI communicator; but
1408  // the function returns an array indexed starting at zero, so we need to
1409  // collect information for each subdomain index anyway, not just for the
1410  // used one.)
1411  const unsigned int n_subdomains =
1412  (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1413  &dof_handler.get_triangulation()) == nullptr ?
1414  [&dof_handler]() {
1415  unsigned int max_subdomain_id = 0;
1416  for (const auto &cell : dof_handler.active_cell_iterators())
1417  max_subdomain_id =
1418  std::max(max_subdomain_id, cell->subdomain_id());
1419  return max_subdomain_id + 1;
1420  }() :
1422  dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1423  &dof_handler.get_triangulation())
1424  ->get_communicator()));
1425  Assert(n_subdomains > *std::max_element(subdomain_association.begin(),
1426  subdomain_association.end()),
1427  ExcInternalError());
1428 
1429  std::vector<::IndexSet> index_sets(
1430  n_subdomains, ::IndexSet(dof_handler.n_dofs()));
1431 
1432  // loop over subdomain_association and populate IndexSet when a
1433  // change in subdomain ID is found
1434  ::types::global_dof_index i_min = 0;
1435  ::types::global_dof_index this_subdomain = subdomain_association[0];
1436 
1437  for (::types::global_dof_index index = 1;
1438  index < subdomain_association.size();
1439  ++index)
1440  {
1441  // found index different from the current one
1442  if (subdomain_association[index] != this_subdomain)
1443  {
1444  index_sets[this_subdomain].add_range(i_min, index);
1445  i_min = index;
1446  this_subdomain = subdomain_association[index];
1447  }
1448  }
1449 
1450  // the very last element is of different index
1451  if (i_min == subdomain_association.size() - 1)
1452  {
1453  index_sets[this_subdomain].add_index(i_min);
1454  }
1455 
1456  // otherwise there are at least two different indices
1457  else
1458  {
1459  index_sets[this_subdomain].add_range(i_min,
1460  subdomain_association.size());
1461  }
1462 
1463  for (unsigned int i = 0; i < n_subdomains; ++i)
1464  index_sets[i].compress();
1465 
1466  return index_sets;
1467  }
1468 
1469  template <int dim, int spacedim>
1470  std::vector<IndexSet>
1472  const DoFHandler<dim, spacedim> &dof_handler)
1473  {
1474  // If the Triangulation is distributed, the only thing we can usefully
1475  // ask is for its locally owned subdomain
1476  Assert((dynamic_cast<
1478  &dof_handler.get_triangulation()) == nullptr),
1479  ExcMessage(
1480  "For parallel::distributed::Triangulation objects and "
1481  "associated DoF handler objects, asking for any information "
1482  "related to a subdomain other than the locally owned one does "
1483  "not make sense."));
1484 
1485  // Collect all the locally owned DoFs
1486  // Note: Even though the distribution of DoFs by the
1487  // locally_owned_dofs_per_subdomain function is pseudo-random, we will
1488  // collect all the DoFs on the subdomain and its layer cell. Therefore, the
1489  // random nature of this function does not play a role in the extraction of
1490  // the locally relevant DoFs
1491  std::vector<IndexSet> dof_set =
1492  locally_owned_dofs_per_subdomain(dof_handler);
1493  const ::types::subdomain_id n_subdomains = dof_set.size();
1494 
1495  // Add the DoFs on the adjacent (equivalent ghost) cells to the IndexSet,
1496  // cache them in a set. Need to check each DoF manually because we can't
1497  // be sure that the DoF range of locally_owned_dofs is really contiguous.
1498  for (::types::subdomain_id subdomain_id = 0;
1499  subdomain_id < n_subdomains;
1500  ++subdomain_id)
1501  {
1502  // Extract the layer of cells around this subdomain
1503  std::function<bool(
1506  const std::vector<
1508  active_halo_layer =
1509  GridTools::compute_active_cell_halo_layer(dof_handler, predicate);
1510 
1511  // Extract DoFs associated with halo layer
1512  std::vector<types::global_dof_index> local_dof_indices;
1513  std::set<types::global_dof_index> subdomain_halo_global_dof_indices;
1514  for (typename std::vector<
1516  const_iterator it_cell = active_halo_layer.begin();
1517  it_cell != active_halo_layer.end();
1518  ++it_cell)
1519  {
1521  &cell = *it_cell;
1522  Assert(
1523  cell->subdomain_id() != subdomain_id,
1524  ExcMessage(
1525  "The subdomain ID of the halo cell should not match that of the vector entry."));
1526 
1527  local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1528  cell->get_dof_indices(local_dof_indices);
1529 
1530  for (const types::global_dof_index local_dof_index :
1531  local_dof_indices)
1532  subdomain_halo_global_dof_indices.insert(local_dof_index);
1533  }
1534 
1535  dof_set[subdomain_id].add_indices(
1536  subdomain_halo_global_dof_indices.begin(),
1537  subdomain_halo_global_dof_indices.end());
1538 
1539  dof_set[subdomain_id].compress();
1540  }
1541 
1542  return dof_set;
1543  }
1544 
1545  template <int dim, int spacedim>
1546  void
1548  const DoFHandler<dim, spacedim> & dof_handler,
1549  std::vector<types::subdomain_id> &subdomain_association)
1550  {
1551  // if the Triangulation is distributed, the only thing we can usefully
1552  // ask is for its locally owned subdomain
1553  Assert((dynamic_cast<
1555  &dof_handler.get_triangulation()) == nullptr),
1556  ExcMessage(
1557  "For parallel::distributed::Triangulation objects and "
1558  "associated DoF handler objects, asking for any subdomain other "
1559  "than the locally owned one does not make sense."));
1560 
1561  Assert(subdomain_association.size() == dof_handler.n_dofs(),
1562  ExcDimensionMismatch(subdomain_association.size(),
1563  dof_handler.n_dofs()));
1564 
1565  // catch an error that happened in some versions of the shared tria
1566  // distribute_dofs() function where we were trying to call this
1567  // function at a point in time when not all internal DoFHandler
1568  // structures were quite set up yet.
1569  Assert(dof_handler.n_dofs() > 0, ExcInternalError());
1570 
1571  // In case this function is executed with parallel::shared::Triangulation
1572  // with possibly artificial cells, we need to take "true" subdomain IDs
1573  // (i.e. without artificial cells). Otherwise we are good to use
1574  // subdomain_id as stored in cell->subdomain_id().
1575  std::vector<types::subdomain_id> cell_owners(
1576  dof_handler.get_triangulation().n_active_cells());
1578  (dynamic_cast<const parallel::shared::Triangulation<dim, spacedim> *>(
1579  &dof_handler.get_triangulation())))
1580  {
1581  cell_owners = tr->get_true_subdomain_ids_of_cells();
1582  Assert(tr->get_true_subdomain_ids_of_cells().size() ==
1583  tr->n_active_cells(),
1584  ExcInternalError());
1585  }
1586  else
1587  {
1588  for (const auto &cell : dof_handler.active_cell_iterators() |
1590  cell_owners[cell->active_cell_index()] = cell->subdomain_id();
1591  }
1592 
1593  // preset all values by an invalid value
1594  std::fill_n(subdomain_association.begin(),
1595  dof_handler.n_dofs(),
1597 
1598  std::vector<types::global_dof_index> local_dof_indices;
1599  local_dof_indices.reserve(
1600  dof_handler.get_fe_collection().max_dofs_per_cell());
1601 
1602  // loop over all cells and record which subdomain a DoF belongs to.
1603  // give to the smaller subdomain_id in case it is on an interface
1604  for (const auto &cell : dof_handler.active_cell_iterators())
1605  {
1607  cell_owners[cell->active_cell_index()];
1608  const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
1609  local_dof_indices.resize(dofs_per_cell);
1610  cell->get_dof_indices(local_dof_indices);
1611 
1612  // set subdomain ids. if dofs already have their values set then
1613  // they must be on partition interfaces. in that case assign them
1614  // to either the previous association or the current processor
1615  // with the smaller subdomain id.
1616  for (unsigned int i = 0; i < dofs_per_cell; ++i)
1617  if (subdomain_association[local_dof_indices[i]] ==
1619  subdomain_association[local_dof_indices[i]] = subdomain_id;
1620  else if (subdomain_association[local_dof_indices[i]] > subdomain_id)
1621  {
1622  subdomain_association[local_dof_indices[i]] = subdomain_id;
1623  }
1624  }
1625 
1626  Assert(std::find(subdomain_association.begin(),
1627  subdomain_association.end(),
1629  subdomain_association.end(),
1630  ExcInternalError());
1631  }
1632 
1633 
1634 
1635  template <int dim, int spacedim>
1636  unsigned int
1638  const DoFHandler<dim, spacedim> &dof_handler,
1639  const types::subdomain_id subdomain)
1640  {
1641  std::vector<types::subdomain_id> subdomain_association(
1642  dof_handler.n_dofs());
1643  get_subdomain_association(dof_handler, subdomain_association);
1644 
1645  return std::count(subdomain_association.begin(),
1646  subdomain_association.end(),
1647  subdomain);
1648  }
1649 
1650 
1651 
1652  template <int dim, int spacedim>
1653  IndexSet
1655  const DoFHandler<dim, spacedim> &dof_handler,
1656  const types::subdomain_id subdomain)
1657  {
1658  // If we have a distributed::Triangulation only allow locally_owned
1659  // subdomain.
1662  (subdomain ==
1663  dof_handler.get_triangulation().locally_owned_subdomain()),
1664  ExcMessage(
1665  "For parallel::distributed::Triangulation objects and "
1666  "associated DoF handler objects, asking for any subdomain other "
1667  "than the locally owned one does not make sense."));
1668 
1669  IndexSet index_set(dof_handler.n_dofs());
1670 
1671  std::vector<types::global_dof_index> local_dof_indices;
1672  local_dof_indices.reserve(
1673  dof_handler.get_fe_collection().max_dofs_per_cell());
1674 
1675  // first generate an unsorted list of all indices which we fill from
1676  // the back. could also insert them directly into the IndexSet, but
1677  // that inserts indices in the middle, which is an O(n^2) algorithm and
1678  // hence too expensive. Could also use std::set, but that is in general
1679  // more expensive than a vector
1680  std::vector<types::global_dof_index> subdomain_indices;
1681 
1682  for (const auto &cell : dof_handler.active_cell_iterators())
1683  if ((cell->is_artificial() == false) &&
1684  (cell->subdomain_id() == subdomain))
1685  {
1686  const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
1687  local_dof_indices.resize(dofs_per_cell);
1688  cell->get_dof_indices(local_dof_indices);
1689  subdomain_indices.insert(subdomain_indices.end(),
1690  local_dof_indices.begin(),
1691  local_dof_indices.end());
1692  }
1693  // sort indices and remove duplicates
1694  std::sort(subdomain_indices.begin(), subdomain_indices.end());
1695  subdomain_indices.erase(std::unique(subdomain_indices.begin(),
1696  subdomain_indices.end()),
1697  subdomain_indices.end());
1698 
1699  // insert into IndexSet
1700  index_set.add_indices(subdomain_indices.begin(), subdomain_indices.end());
1701  index_set.compress();
1702 
1703  return index_set;
1704  }
1705 
1706 
1707 
1708  template <int dim, int spacedim>
1709  void
1711  const DoFHandler<dim, spacedim> &dof_handler,
1712  const types::subdomain_id subdomain,
1713  std::vector<unsigned int> & n_dofs_on_subdomain)
1714  {
1715  Assert(n_dofs_on_subdomain.size() == dof_handler.get_fe(0).n_components(),
1716  ExcDimensionMismatch(n_dofs_on_subdomain.size(),
1717  dof_handler.get_fe(0).n_components()));
1718  std::fill(n_dofs_on_subdomain.begin(), n_dofs_on_subdomain.end(), 0);
1719 
1720  // Make sure there are at least some cells with this subdomain id
1721  Assert(std::any_of(
1722  dof_handler.begin_active(),
1724  dof_handler.end()},
1725  [subdomain](
1726  const typename DoFHandler<dim, spacedim>::cell_accessor &cell) {
1727  return cell.subdomain_id() == subdomain;
1728  }),
1729  ExcMessage("There are no cells for the given subdomain!"));
1730 
1731  std::vector<types::subdomain_id> subdomain_association(
1732  dof_handler.n_dofs());
1733  get_subdomain_association(dof_handler, subdomain_association);
1734 
1735  std::vector<unsigned char> component_association(dof_handler.n_dofs());
1737  std::vector<bool>(),
1738  component_association);
1739 
1740  for (unsigned int c = 0; c < dof_handler.get_fe(0).n_components(); ++c)
1741  {
1742  for (types::global_dof_index i = 0; i < dof_handler.n_dofs(); ++i)
1743  if ((subdomain_association[i] == subdomain) &&
1744  (component_association[i] == static_cast<unsigned char>(c)))
1745  ++n_dofs_on_subdomain[c];
1746  }
1747  }
1748 
1749 
1750 
1751  namespace internal
1752  {
1753  // TODO: why is this function so complicated? It would be nice to have
1754  // comments that explain why we can't just loop over all components and
1755  // count the entries in dofs_by_component that have this component's
1756  // index
1757  template <int dim, int spacedim>
1758  void
1760  const std::vector<unsigned char> & dofs_by_component,
1761  const std::vector<unsigned int> & target_component,
1762  const bool only_once,
1763  std::vector<types::global_dof_index> &dofs_per_component,
1764  unsigned int & component)
1765  {
1766  for (unsigned int b = 0; b < fe.n_base_elements(); ++b)
1767  {
1768  const FiniteElement<dim, spacedim> &base = fe.base_element(b);
1769  // Dimension of base element
1770  unsigned int d = base.n_components();
1771 
1772  for (unsigned int m = 0; m < fe.element_multiplicity(b); ++m)
1773  {
1774  if (base.n_base_elements() > 1)
1775  resolve_components(base,
1776  dofs_by_component,
1777  target_component,
1778  only_once,
1779  dofs_per_component,
1780  component);
1781  else
1782  {
1783  for (unsigned int dd = 0; dd < d; ++dd, ++component)
1784  dofs_per_component[target_component[component]] +=
1785  std::count(dofs_by_component.begin(),
1786  dofs_by_component.end(),
1787  component);
1788 
1789  // if we have non-primitive FEs and want all components
1790  // to show the number of dofs, need to copy the result to
1791  // those components
1792  if (!base.is_primitive() && !only_once)
1793  for (unsigned int dd = 1; dd < d; ++dd)
1794  dofs_per_component[target_component[component - d + dd]] =
1795  dofs_per_component[target_component[component - d]];
1796  }
1797  }
1798  }
1799  }
1800 
1801 
1802  template <int dim, int spacedim>
1803  void
1805  const std::vector<unsigned char> & dofs_by_component,
1806  const std::vector<unsigned int> & target_component,
1807  const bool only_once,
1808  std::vector<types::global_dof_index> &dofs_per_component,
1809  unsigned int & component)
1810  {
1811  // assert that all elements in the collection have the same structure
1812  // (base elements and multiplicity, components per base element) and
1813  // then simply call the function above
1814  for (unsigned int fe = 1; fe < fe_collection.size(); ++fe)
1815  {
1816  Assert(fe_collection[fe].n_components() ==
1817  fe_collection[0].n_components(),
1818  ExcNotImplemented());
1819  Assert(fe_collection[fe].n_base_elements() ==
1820  fe_collection[0].n_base_elements(),
1821  ExcNotImplemented());
1822  for (unsigned int b = 0; b < fe_collection[0].n_base_elements(); ++b)
1823  {
1824  Assert(fe_collection[fe].base_element(b).n_components() ==
1825  fe_collection[0].base_element(b).n_components(),
1826  ExcNotImplemented());
1827  Assert(fe_collection[fe].base_element(b).n_base_elements() ==
1828  fe_collection[0].base_element(b).n_base_elements(),
1829  ExcNotImplemented());
1830  }
1831  }
1832 
1833  resolve_components(fe_collection[0],
1834  dofs_by_component,
1835  target_component,
1836  only_once,
1837  dofs_per_component,
1838  component);
1839  }
1840  } // namespace internal
1841 
1842 
1843 
1844  namespace internal
1845  {
1846  namespace
1847  {
1851  template <int dim, int spacedim>
1852  bool
1853  all_elements_are_primitive(const FiniteElement<dim, spacedim> &fe)
1854  {
1855  return fe.is_primitive();
1856  }
1857 
1858 
1863  template <int dim, int spacedim>
1864  bool
1865  all_elements_are_primitive(
1866  const ::hp::FECollection<dim, spacedim> &fe_collection)
1867  {
1868  for (unsigned int i = 0; i < fe_collection.size(); ++i)
1869  if (fe_collection[i].is_primitive() == false)
1870  return false;
1871 
1872  return true;
1873  }
1874  } // namespace
1875  } // namespace internal
1876 
1877 
1878 
1879  template <int dim, int spacedim>
1880  std::vector<types::global_dof_index>
1882  const DoFHandler<dim, spacedim> &dof_handler,
1883  const bool only_once,
1884  const std::vector<unsigned int> &target_component_)
1885  {
1886  const unsigned int n_components = dof_handler.get_fe(0).n_components();
1887 
1888  // If the empty vector was given as default argument, set up this
1889  // vector as identity.
1890  std::vector<unsigned int> target_component = target_component_;
1891  if (target_component.size() == 0)
1892  {
1893  target_component.resize(n_components);
1894  for (unsigned int i = 0; i < n_components; ++i)
1895  target_component[i] = i;
1896  }
1897  else
1898  Assert(target_component.size() == n_components,
1899  ExcDimensionMismatch(target_component.size(), n_components));
1900 
1901 
1902  const unsigned int max_component =
1903  *std::max_element(target_component.begin(), target_component.end());
1904  const unsigned int n_target_components = max_component + 1;
1905 
1906  std::vector<types::global_dof_index> dofs_per_component(
1907  n_target_components, types::global_dof_index(0));
1908 
1909  // special case for only one component. treat this first since it does
1910  // not require any computations
1911  if (n_components == 1)
1912  {
1913  dofs_per_component[0] = dof_handler.n_locally_owned_dofs();
1914  return dofs_per_component;
1915  }
1916 
1917 
1918  // otherwise determine the number of dofs in each component separately.
1919  // do so in parallel
1920  std::vector<unsigned char> dofs_by_component(
1921  dof_handler.n_locally_owned_dofs());
1923  ComponentMask(),
1924  dofs_by_component);
1925 
1926  // next count what we got
1927  unsigned int component = 0;
1929  dofs_by_component,
1930  target_component,
1931  only_once,
1932  dofs_per_component,
1933  component);
1934  Assert(n_components == component, ExcInternalError());
1935 
1936  // finally sanity check. this is only valid if the finite element is
1937  // actually primitive, so exclude other elements from this
1938  Assert((internal::all_elements_are_primitive(
1939  dof_handler.get_fe_collection()) == false) ||
1940  (std::accumulate(dofs_per_component.begin(),
1941  dofs_per_component.end(),
1943  dof_handler.n_locally_owned_dofs()),
1944  ExcInternalError());
1945 
1946  // reduce information from all CPUs
1947 #ifdef DEAL_II_WITH_MPI
1948 
1950  (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1951  &dof_handler.get_triangulation())))
1952  {
1953  std::vector<types::global_dof_index> local_dof_count =
1954  dofs_per_component;
1955 
1956  const int ierr = MPI_Allreduce(local_dof_count.data(),
1957  dofs_per_component.data(),
1958  n_target_components,
1960  MPI_SUM,
1961  tria->get_communicator());
1962  AssertThrowMPI(ierr);
1963  }
1964 #endif
1965 
1966  return dofs_per_component;
1967  }
1968 
1969 
1970 
1971  template <int dim, int spacedim>
1972  std::vector<types::global_dof_index>
1974  const std::vector<unsigned int> &target_block_)
1975  {
1976  const ::hp::FECollection<dim, spacedim> &fe_collection =
1977  dof_handler.get_fe_collection();
1978  Assert(fe_collection.size() < 256, ExcNotImplemented());
1979 
1980  // If the empty vector for target_block(e.g., as default argument), then
1981  // set up this vector as identity. We do this set up with the first
1982  // element of the collection, but the whole thing can only work if
1983  // all elements have the same number of blocks anyway -- so check
1984  // that right after
1985  const unsigned int n_blocks = fe_collection[0].n_blocks();
1986 
1987  std::vector<unsigned int> target_block = target_block_;
1988  if (target_block.size() == 0)
1989  {
1990  target_block.resize(fe_collection[0].n_blocks());
1991  for (unsigned int i = 0; i < n_blocks; ++i)
1992  target_block[i] = i;
1993  }
1994  else
1995  Assert(target_block.size() == n_blocks,
1996  ExcDimensionMismatch(target_block.size(), n_blocks));
1997  for (unsigned int f = 1; f < fe_collection.size(); ++f)
1998  Assert(fe_collection[0].n_blocks() == fe_collection[f].n_blocks(),
1999  ExcMessage("This function can only work if all elements in a "
2000  "collection have the same number of blocks."));
2001 
2002  // special case for only one block. treat this first since it does
2003  // not require any computations
2004  if (n_blocks == 1)
2005  {
2006  std::vector<types::global_dof_index> dofs_per_block(1);
2007  dofs_per_block[0] = dof_handler.n_dofs();
2008  return dofs_per_block;
2009  }
2010 
2011  // Otherwise set up the right-sized object and start working
2012  const unsigned int max_block =
2013  *std::max_element(target_block.begin(), target_block.end());
2014  const unsigned int n_target_blocks = max_block + 1;
2015 
2016  std::vector<types::global_dof_index> dofs_per_block(n_target_blocks);
2017 
2018  // Loop over the elements of the collection, but really only consider
2019  // the last element (see #9271)
2020  for (unsigned int this_fe = fe_collection.size() - 1;
2021  this_fe < fe_collection.size();
2022  ++this_fe)
2023  {
2024  const FiniteElement<dim, spacedim> &fe = fe_collection[this_fe];
2025 
2026  std::vector<unsigned char> dofs_by_block(
2027  dof_handler.n_locally_owned_dofs());
2028  internal::get_block_association(dof_handler, dofs_by_block);
2029 
2030  // next count what we got
2031  for (unsigned int block = 0; block < fe.n_blocks(); ++block)
2032  dofs_per_block[target_block[block]] +=
2033  std::count(dofs_by_block.begin(), dofs_by_block.end(), block);
2034 
2035 #ifdef DEAL_II_WITH_MPI
2036  // if we are working on a parallel mesh, we now need to collect
2037  // this information from all processors
2039  (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
2040  &dof_handler.get_triangulation())))
2041  {
2042  std::vector<types::global_dof_index> local_dof_count =
2043  dofs_per_block;
2044  const int ierr = MPI_Allreduce(local_dof_count.data(),
2045  dofs_per_block.data(),
2046  n_target_blocks,
2048  MPI_SUM,
2049  tria->get_communicator());
2050  AssertThrowMPI(ierr);
2051  }
2052 #endif
2053  }
2054 
2055  return dofs_per_block;
2056  }
2057 
2058 
2059 
2060  template <int dim, int spacedim>
2061  void
2063  std::vector<types::global_dof_index> &mapping)
2064  {
2065  mapping.clear();
2066  mapping.insert(mapping.end(),
2067  dof_handler.n_dofs(),
2069 
2070  std::vector<types::global_dof_index> dofs_on_face;
2071  dofs_on_face.reserve(dof_handler.get_fe_collection().max_dofs_per_face());
2072  types::global_dof_index next_boundary_index = 0;
2073 
2074  // now loop over all cells and check whether their faces are at the
2075  // boundary. note that we need not take special care of single lines
2076  // being at the boundary (using @p{cell->has_boundary_lines}), since we
2077  // do not support boundaries of dimension dim-2, and so every isolated
2078  // boundary line is also part of a boundary face which we will be
2079  // visiting sooner or later
2080  for (const auto &cell : dof_handler.active_cell_iterators())
2081  for (const unsigned int f : cell->face_indices())
2082  if (cell->at_boundary(f))
2083  {
2084  const unsigned int dofs_per_face =
2085  cell->get_fe().n_dofs_per_face(f);
2086  dofs_on_face.resize(dofs_per_face);
2087  cell->face(f)->get_dof_indices(dofs_on_face,
2088  cell->active_fe_index());
2089  for (unsigned int i = 0; i < dofs_per_face; ++i)
2090  if (mapping[dofs_on_face[i]] == numbers::invalid_dof_index)
2091  mapping[dofs_on_face[i]] = next_boundary_index++;
2092  }
2093 
2094  AssertDimension(next_boundary_index, dof_handler.n_boundary_dofs());
2095  }
2096 
2097 
2098 
2099  template <int dim, int spacedim>
2100  void
2102  const std::set<types::boundary_id> &boundary_ids,
2103  std::vector<types::global_dof_index> &mapping)
2104  {
2105  Assert(boundary_ids.find(numbers::internal_face_boundary_id) ==
2106  boundary_ids.end(),
2108 
2109  mapping.clear();
2110  mapping.insert(mapping.end(),
2111  dof_handler.n_dofs(),
2113 
2114  // return if there is nothing to do
2115  if (boundary_ids.size() == 0)
2116  return;
2117 
2118  std::vector<types::global_dof_index> dofs_on_face;
2119  dofs_on_face.reserve(dof_handler.get_fe_collection().max_dofs_per_face());
2120  types::global_dof_index next_boundary_index = 0;
2121 
2122  for (const auto &cell : dof_handler.active_cell_iterators())
2123  for (const unsigned int f : cell->face_indices())
2124  if (boundary_ids.find(cell->face(f)->boundary_id()) !=
2125  boundary_ids.end())
2126  {
2127  const unsigned int dofs_per_face =
2128  cell->get_fe().n_dofs_per_face(f);
2129  dofs_on_face.resize(dofs_per_face);
2130  cell->face(f)->get_dof_indices(dofs_on_face,
2131  cell->active_fe_index());
2132  for (unsigned int i = 0; i < dofs_per_face; ++i)
2133  if (mapping[dofs_on_face[i]] == numbers::invalid_dof_index)
2134  mapping[dofs_on_face[i]] = next_boundary_index++;
2135  }
2136 
2137  AssertDimension(next_boundary_index,
2138  dof_handler.n_boundary_dofs(boundary_ids));
2139  }
2140 
2141  namespace internal
2142  {
2143  namespace
2144  {
2145  template <int dim, int spacedim>
2146  void
2148  const hp::MappingCollection<dim, spacedim> & mapping,
2149  const DoFHandler<dim, spacedim> & dof_handler,
2150  std::map<types::global_dof_index, Point<spacedim>> &support_points,
2151  const ComponentMask & in_mask)
2152  {
2153  const hp::FECollection<dim, spacedim> &fe_collection =
2154  dof_handler.get_fe_collection();
2155  hp::QCollection<dim> q_coll_dummy;
2156 
2157  for (unsigned int fe_index = 0; fe_index < fe_collection.size();
2158  ++fe_index)
2159  {
2160  // check whether every FE in the collection has support points
2161  Assert((fe_collection[fe_index].n_dofs_per_cell() == 0) ||
2162  (fe_collection[fe_index].has_support_points()),
2164  q_coll_dummy.push_back(Quadrature<dim>(
2165  fe_collection[fe_index].get_unit_support_points()));
2166  }
2167 
2168  // Take care of components
2169  const ComponentMask mask =
2170  (in_mask.size() == 0 ?
2171  ComponentMask(fe_collection.n_components(), true) :
2172  in_mask);
2173 
2174  // Now loop over all cells and enquire the support points on each
2175  // of these. we use dummy quadrature formulas where the quadrature
2176  // points are located at the unit support points to enquire the
2177  // location of the support points in real space.
2178  //
2179  // The weights of the quadrature rule have been set to invalid
2180  // values by the used constructor.
2181  hp::FEValues<dim, spacedim> hp_fe_values(mapping,
2182  fe_collection,
2183  q_coll_dummy,
2185 
2186  std::vector<types::global_dof_index> local_dof_indices;
2187  for (const auto &cell : dof_handler.active_cell_iterators())
2188  // only work on locally relevant cells
2189  if (cell->is_artificial() == false)
2190  {
2191  hp_fe_values.reinit(cell);
2192  const FEValues<dim, spacedim> &fe_values =
2193  hp_fe_values.get_present_fe_values();
2194 
2195  local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
2196  cell->get_dof_indices(local_dof_indices);
2197 
2198  const std::vector<Point<spacedim>> &points =
2199  fe_values.get_quadrature_points();
2200  for (unsigned int i = 0; i < cell->get_fe().n_dofs_per_cell();
2201  ++i)
2202  {
2203  const unsigned int dof_comp =
2204  cell->get_fe().system_to_component_index(i).first;
2205 
2206  // insert the values into the map if it is a valid component
2207  if (mask[dof_comp])
2208  support_points[local_dof_indices[i]] = points[i];
2209  }
2210  }
2211  }
2212 
2213 
2214  template <int dim, int spacedim>
2215  void
2217  const hp::MappingCollection<dim, spacedim> &mapping,
2218  const DoFHandler<dim, spacedim> & dof_handler,
2219  std::vector<Point<spacedim>> & support_points,
2220  const ComponentMask & mask)
2221  {
2222  // get the data in the form of the map as above
2223  std::map<types::global_dof_index, Point<spacedim>> x_support_points;
2225  dof_handler,
2226  x_support_points,
2227  mask);
2228 
2229  // now convert from the map to the linear vector. make sure every
2230  // entry really appeared in the map
2231  for (types::global_dof_index i = 0; i < dof_handler.n_dofs(); ++i)
2232  {
2233  Assert(x_support_points.find(i) != x_support_points.end(),
2234  ExcInternalError());
2235 
2236  support_points[i] = x_support_points[i];
2237  }
2238  }
2239  } // namespace
2240  } // namespace internal
2241 
2242  template <int dim, int spacedim>
2243  void
2245  const DoFHandler<dim, spacedim> &dof_handler,
2246  std::vector<Point<spacedim>> & support_points,
2247  const ComponentMask & mask)
2248  {
2249  AssertDimension(support_points.size(), dof_handler.n_dofs());
2250  Assert((dynamic_cast<
2252  &dof_handler.get_triangulation()) == nullptr),
2253  ExcMessage(
2254  "This function can not be used with distributed triangulations. "
2255  "See the documentation for more information."));
2256 
2257  // Let the internal function do all the work, just make sure that it
2258  // gets a MappingCollection
2259  const hp::MappingCollection<dim, spacedim> mapping_collection(mapping);
2260 
2261  internal::map_dofs_to_support_points(mapping_collection,
2262  dof_handler,
2263  support_points,
2264  mask);
2265  }
2266 
2267 
2268  template <int dim, int spacedim>
2269  void
2271  const hp::MappingCollection<dim, spacedim> &mapping,
2272  const DoFHandler<dim, spacedim> & dof_handler,
2273  std::vector<Point<spacedim>> & support_points,
2274  const ComponentMask & mask)
2275  {
2276  AssertDimension(support_points.size(), dof_handler.n_dofs());
2277  Assert((dynamic_cast<
2279  &dof_handler.get_triangulation()) == nullptr),
2280  ExcMessage(
2281  "This function can not be used with distributed triangulations. "
2282  "See the documentation for more information."));
2283 
2284  // Let the internal function do all the work, just make sure that it
2285  // gets a MappingCollection
2287  dof_handler,
2288  support_points,
2289  mask);
2290  }
2291 
2292 
2293  template <int dim, int spacedim>
2294  void
2296  const Mapping<dim, spacedim> & mapping,
2297  const DoFHandler<dim, spacedim> & dof_handler,
2298  std::map<types::global_dof_index, Point<spacedim>> &support_points,
2299  const ComponentMask & mask)
2300  {
2301  support_points.clear();
2302 
2303  // Let the internal function do all the work, just make sure that it
2304  // gets a MappingCollection
2305  const hp::MappingCollection<dim, spacedim> mapping_collection(mapping);
2306 
2307  internal::map_dofs_to_support_points(mapping_collection,
2308  dof_handler,
2309  support_points,
2310  mask);
2311  }
2312 
2313 
2314  template <int dim, int spacedim>
2315  void
2317  const hp::MappingCollection<dim, spacedim> & mapping,
2318  const DoFHandler<dim, spacedim> & dof_handler,
2319  std::map<types::global_dof_index, Point<spacedim>> &support_points,
2320  const ComponentMask & mask)
2321  {
2322  support_points.clear();
2323 
2324  // Let the internal function do all the work, just make sure that it
2325  // gets a MappingCollection
2327  dof_handler,
2328  support_points,
2329  mask);
2330  }
2331 
2332  template <int spacedim>
2333  void
2335  std::ostream & out,
2336  const std::map<types::global_dof_index, Point<spacedim>> &support_points)
2337  {
2338  AssertThrow(out.fail() == false, ExcIO());
2339 
2340  std::map<Point<spacedim>,
2341  std::vector<types::global_dof_index>,
2343  point_map;
2344 
2345  // convert to map point -> list of DoFs
2346  for (const auto &it : support_points)
2347  {
2348  std::vector<types::global_dof_index> &v = point_map[it.second];
2349  v.push_back(it.first);
2350  }
2351 
2352  // print the newly created map:
2353  for (const auto &it : point_map)
2354  {
2355  out << it.first << " \"";
2356  const std::vector<types::global_dof_index> &v = it.second;
2357  for (unsigned int i = 0; i < v.size(); ++i)
2358  {
2359  if (i > 0)
2360  out << ", ";
2361  out << v[i];
2362  }
2363 
2364  out << "\"\n";
2365  }
2366 
2367  out << std::flush;
2368  }
2369 
2370 
2371  template <int dim, int spacedim>
2372  void
2374  const Table<2, Coupling> & table,
2375  std::vector<Table<2, Coupling>> &tables_by_block)
2376  {
2377  if (dof_handler.has_hp_capabilities() == false)
2378  {
2379  const FiniteElement<dim, spacedim> &fe = dof_handler.get_fe();
2380  const unsigned int nb = fe.n_blocks();
2381 
2382  tables_by_block.resize(1);
2383  tables_by_block[0].reinit(nb, nb);
2384  tables_by_block[0].fill(none);
2385 
2386  for (unsigned int i = 0; i < fe.n_components(); ++i)
2387  {
2388  const unsigned int ib = fe.component_to_block_index(i);
2389  for (unsigned int j = 0; j < fe.n_components(); ++j)
2390  {
2391  const unsigned int jb = fe.component_to_block_index(j);
2392  tables_by_block[0](ib, jb) |= table(i, j);
2393  }
2394  }
2395  }
2396  else
2397  {
2398  const hp::FECollection<dim> &fe_collection =
2399  dof_handler.get_fe_collection();
2400  tables_by_block.resize(fe_collection.size());
2401 
2402  for (unsigned int f = 0; f < fe_collection.size(); ++f)
2403  {
2404  const FiniteElement<dim, spacedim> &fe = fe_collection[f];
2405 
2406  const unsigned int nb = fe.n_blocks();
2407  tables_by_block[f].reinit(nb, nb);
2408  tables_by_block[f].fill(none);
2409  for (unsigned int i = 0; i < fe.n_components(); ++i)
2410  {
2411  const unsigned int ib = fe.component_to_block_index(i);
2412  for (unsigned int j = 0; j < fe.n_components(); ++j)
2413  {
2414  const unsigned int jb = fe.component_to_block_index(j);
2415  tables_by_block[f](ib, jb) |= table(i, j);
2416  }
2417  }
2418  }
2419  }
2420  }
2421 
2422 
2423 
2424  template <int dim, int spacedim>
2425  void
2427  const DoFHandler<dim, spacedim> &dof_handler,
2428  const unsigned int level,
2429  const std::vector<bool> & selected_dofs,
2430  const types::global_dof_index offset)
2431  {
2432  std::vector<types::global_dof_index> indices;
2433 
2434  unsigned int i = 0;
2435 
2436  for (const auto &cell : dof_handler.cell_iterators_on_level(level))
2437  if (cell->is_locally_owned_on_level())
2438  ++i;
2439  block_list.reinit(i,
2440  dof_handler.n_dofs(),
2441  dof_handler.get_fe().n_dofs_per_cell());
2442  i = 0;
2443  for (const auto &cell : dof_handler.cell_iterators_on_level(level))
2444  if (cell->is_locally_owned_on_level())
2445  {
2446  indices.resize(cell->get_fe().n_dofs_per_cell());
2447  cell->get_mg_dof_indices(indices);
2448 
2449  if (selected_dofs.size() != 0)
2450  AssertDimension(indices.size(), selected_dofs.size());
2451 
2452  for (types::global_dof_index j = 0; j < indices.size(); ++j)
2453  {
2454  if (selected_dofs.size() == 0)
2455  block_list.add(i, indices[j] - offset);
2456  else
2457  {
2458  if (selected_dofs[j])
2459  block_list.add(i, indices[j] - offset);
2460  }
2461  }
2462  ++i;
2463  }
2464  }
2465 
2466 
2467  template <int dim, int spacedim>
2468  void
2470  const DoFHandler<dim, spacedim> &dof_handler,
2471  const unsigned int level,
2472  const bool interior_only)
2473  {
2474  const FiniteElement<dim> &fe = dof_handler.get_fe();
2475  block_list.reinit(1, dof_handler.n_dofs(level), dof_handler.n_dofs(level));
2476 
2477  std::vector<types::global_dof_index> indices;
2478  std::vector<bool> exclude;
2479 
2480  for (const auto &cell : dof_handler.cell_iterators_on_level(level))
2481  {
2482  indices.resize(cell->get_fe().n_dofs_per_cell());
2483  cell->get_mg_dof_indices(indices);
2484 
2485  if (interior_only)
2486  {
2487  // Exclude degrees of freedom on faces opposite to the vertex
2488  exclude.resize(fe.n_dofs_per_cell());
2489  std::fill(exclude.begin(), exclude.end(), false);
2490 
2491  for (const unsigned int face : cell->face_indices())
2492  if (cell->at_boundary(face) ||
2493  cell->neighbor(face)->level() != cell->level())
2494  for (unsigned int i = 0; i < fe.n_dofs_per_face(face); ++i)
2495  exclude[fe.face_to_cell_index(i, face)] = true;
2496  for (types::global_dof_index j = 0; j < indices.size(); ++j)
2497  if (!exclude[j])
2498  block_list.add(0, indices[j]);
2499  }
2500  else
2501  {
2502  for (const auto index : indices)
2503  block_list.add(0, index);
2504  }
2505  }
2506  }
2507 
2508 
2509  template <int dim, int spacedim>
2510  void
2512  const DoFHandler<dim, spacedim> &dof_handler,
2513  const unsigned int level,
2514  const bool interior_dofs_only,
2515  const bool boundary_dofs)
2516  {
2517  Assert(level > 0 && level < dof_handler.get_triangulation().n_levels(),
2518  ExcIndexRange(level, 1, dof_handler.get_triangulation().n_levels()));
2519 
2520  std::vector<types::global_dof_index> indices;
2521  std::vector<bool> exclude;
2522 
2523  unsigned int block = 0;
2524  for (const auto &pcell : dof_handler.cell_iterators_on_level(level - 1))
2525  {
2526  if (pcell->is_active())
2527  continue;
2528 
2529  for (unsigned int child = 0; child < pcell->n_children(); ++child)
2530  {
2531  const auto cell = pcell->child(child);
2532 
2533  // For hp, only this line here would have to be replaced.
2534  const FiniteElement<dim> &fe = dof_handler.get_fe();
2535  const unsigned int n_dofs = fe.n_dofs_per_cell();
2536  indices.resize(n_dofs);
2537  exclude.resize(n_dofs);
2538  std::fill(exclude.begin(), exclude.end(), false);
2539  cell->get_mg_dof_indices(indices);
2540 
2541  if (interior_dofs_only)
2542  {
2543  // Eliminate dofs on faces of the child which are on faces
2544  // of the parent
2545  for (unsigned int d = 0; d < dim; ++d)
2546  {
2547  const unsigned int face =
2549  for (unsigned int i = 0; i < fe.n_dofs_per_face(face); ++i)
2550  exclude[fe.face_to_cell_index(i, face)] = true;
2551  }
2552 
2553  // Now remove all degrees of freedom on the domain boundary
2554  // from the exclusion list
2555  if (boundary_dofs)
2556  for (const unsigned int face :
2558  if (cell->at_boundary(face))
2559  for (unsigned int i = 0; i < fe.n_dofs_per_face(face);
2560  ++i)
2561  exclude[fe.face_to_cell_index(i, face)] = false;
2562  }
2563 
2564  for (unsigned int i = 0; i < n_dofs; ++i)
2565  if (!exclude[i])
2566  block_list.add(block, indices[i]);
2567  }
2568  ++block;
2569  }
2570  }
2571 
2572  template <int dim, int spacedim>
2573  std::vector<unsigned int>
2575  const DoFHandler<dim, spacedim> &dof_handler,
2576  const unsigned int level,
2577  const bool interior_only,
2578  const bool boundary_patches,
2579  const bool level_boundary_patches,
2580  const bool single_cell_patches,
2581  const bool invert_vertex_mapping)
2582  {
2583  const unsigned int n_blocks = dof_handler.get_fe().n_blocks();
2584  BlockMask exclude_boundary_dofs = BlockMask(n_blocks, interior_only);
2585  return make_vertex_patches(block_list,
2586  dof_handler,
2587  level,
2588  exclude_boundary_dofs,
2589  boundary_patches,
2590  level_boundary_patches,
2591  single_cell_patches,
2592  invert_vertex_mapping);
2593  }
2594 
2595  template <int dim, int spacedim>
2596  std::vector<unsigned int>
2598  const DoFHandler<dim, spacedim> &dof_handler,
2599  const unsigned int level,
2600  const BlockMask & exclude_boundary_dofs,
2601  const bool boundary_patches,
2602  const bool level_boundary_patches,
2603  const bool single_cell_patches,
2604  const bool invert_vertex_mapping)
2605  {
2606  // Vector mapping from vertex index in the triangulation to consecutive
2607  // block indices on this level The number of cells at a vertex
2608  std::vector<unsigned int> vertex_cell_count(
2609  dof_handler.get_triangulation().n_vertices(), 0);
2610 
2611  // Is a vertex at the boundary?
2612  std::vector<bool> vertex_boundary(
2613  dof_handler.get_triangulation().n_vertices(), false);
2614 
2615  std::vector<unsigned int> vertex_mapping(
2616  dof_handler.get_triangulation().n_vertices(),
2618 
2619  // Estimate for the number of dofs at this point
2620  std::vector<unsigned int> vertex_dof_count(
2621  dof_handler.get_triangulation().n_vertices(), 0);
2622 
2623  // Identify all vertices active on this level and remember some data
2624  // about them
2625  for (const auto &cell : dof_handler.cell_iterators_on_level(level))
2626  for (const unsigned int v : cell->vertex_indices())
2627  {
2628  const unsigned int vg = cell->vertex_index(v);
2629  vertex_dof_count[vg] += cell->get_fe().n_dofs_per_cell();
2630  ++vertex_cell_count[vg];
2631  for (unsigned int d = 0; d < dim; ++d)
2632  {
2633  const unsigned int face = GeometryInfo<dim>::vertex_to_face[v][d];
2634  if (cell->at_boundary(face))
2635  vertex_boundary[vg] = true;
2636  else if ((!level_boundary_patches) &&
2637  (cell->neighbor(face)->level() !=
2638  static_cast<int>(level)))
2639  vertex_boundary[vg] = true;
2640  }
2641  }
2642  // From now on, only vertices with positive dof count are "in".
2643 
2644  // Remove vertices at boundaries or in corners
2645  for (unsigned int vg = 0; vg < vertex_dof_count.size(); ++vg)
2646  if ((!single_cell_patches && vertex_cell_count[vg] < 2) ||
2647  (!boundary_patches && vertex_boundary[vg]))
2648  vertex_dof_count[vg] = 0;
2649 
2650  // Create a mapping from all vertices to the ones used here
2651  unsigned int n_vertex_count = 0;
2652  for (unsigned int vg = 0; vg < vertex_mapping.size(); ++vg)
2653  if (vertex_dof_count[vg] != 0)
2654  vertex_mapping[vg] = n_vertex_count++;
2655 
2656  // Compactify dof count
2657  for (unsigned int vg = 0; vg < vertex_mapping.size(); ++vg)
2658  if (vertex_dof_count[vg] != 0)
2659  vertex_dof_count[vertex_mapping[vg]] = vertex_dof_count[vg];
2660 
2661  // Now that we have all the data, we reduce it to the part we actually
2662  // want
2663  vertex_dof_count.resize(n_vertex_count);
2664 
2665  // At this point, the list of patches is ready. Now we enter the dofs
2666  // into the sparsity pattern.
2667  block_list.reinit(vertex_dof_count.size(),
2668  dof_handler.n_dofs(level),
2669  vertex_dof_count);
2670 
2671  std::vector<types::global_dof_index> indices;
2672  std::vector<bool> exclude;
2673 
2674  for (const auto &cell : dof_handler.cell_iterators_on_level(level))
2675  {
2676  const FiniteElement<dim> &fe = cell->get_fe();
2677  indices.resize(fe.n_dofs_per_cell());
2678  cell->get_mg_dof_indices(indices);
2679 
2680  for (const unsigned int v : cell->vertex_indices())
2681  {
2682  const unsigned int vg = cell->vertex_index(v);
2683  const unsigned int block = vertex_mapping[vg];
2684  if (block == numbers::invalid_unsigned_int)
2685  continue;
2686 
2687  // Collect excluded dofs for some block(s) if boundary dofs
2688  // for a block are decided to be excluded
2689  if (exclude_boundary_dofs.size() == 0 ||
2690  exclude_boundary_dofs.n_selected_blocks() != 0)
2691  {
2692  // Exclude degrees of freedom on faces opposite to the
2693  // vertex
2694  exclude.resize(fe.n_dofs_per_cell());
2695  std::fill(exclude.begin(), exclude.end(), false);
2696 
2697  for (unsigned int d = 0; d < dim; ++d)
2698  {
2699  const unsigned int a_face =
2701  const unsigned int face =
2703  for (unsigned int i = 0; i < fe.n_dofs_per_face(face); ++i)
2704  {
2705  // For each dof, get the block it is in and decide to
2706  // exclude it or not
2707  if (exclude_boundary_dofs[fe.system_to_block_index(
2708  fe.face_to_cell_index(
2709  i, face))
2710  .first] == true)
2711  exclude[fe.face_to_cell_index(i, face)] = true;
2712  }
2713  }
2714  for (unsigned int j = 0; j < indices.size(); ++j)
2715  if (!exclude[j])
2716  block_list.add(block, indices[j]);
2717  }
2718  else
2719  {
2720  for (const auto index : indices)
2721  block_list.add(block, index);
2722  }
2723  }
2724  }
2725 
2726  if (invert_vertex_mapping)
2727  {
2728  // Compress vertex mapping
2729  unsigned int n_vertex_count = 0;
2730  for (unsigned int vg = 0; vg < vertex_mapping.size(); ++vg)
2731  if (vertex_mapping[vg] != numbers::invalid_unsigned_int)
2732  vertex_mapping[n_vertex_count++] = vg;
2733 
2734  // Now we reduce it to the part we actually want
2735  vertex_mapping.resize(n_vertex_count);
2736  }
2737 
2738  return vertex_mapping;
2739  }
2740 
2741 
2742  template <int dim, int spacedim>
2743  unsigned int
2745  const std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
2746  &patch)
2747  {
2748  std::set<types::global_dof_index> dofs_on_patch;
2749  std::vector<types::global_dof_index> local_dof_indices;
2750 
2751  // loop over the cells in the patch and get the DoFs on each.
2752  // add all of them to a std::set which automatically makes sure
2753  // all duplicates are ignored
2754  for (unsigned int i = 0; i < patch.size(); ++i)
2755  {
2757  patch[i];
2758  Assert(cell->is_artificial() == false,
2759  ExcMessage("This function can not be called with cells that are "
2760  "not either locally owned or ghost cells."));
2761  local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
2762  cell->get_dof_indices(local_dof_indices);
2763  dofs_on_patch.insert(local_dof_indices.begin(),
2764  local_dof_indices.end());
2765  }
2766 
2767  // now return the number of DoFs (duplicates were ignored)
2768  return dofs_on_patch.size();
2769  }
2770 
2771 
2772 
2773  template <int dim, int spacedim>
2774  std::vector<types::global_dof_index>
2776  const std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
2777  &patch)
2778  {
2779  std::set<types::global_dof_index> dofs_on_patch;
2780  std::vector<types::global_dof_index> local_dof_indices;
2781 
2782  // loop over the cells in the patch and get the DoFs on each.
2783  // add all of them to a std::set which automatically makes sure
2784  // all duplicates are ignored
2785  for (unsigned int i = 0; i < patch.size(); ++i)
2786  {
2788  patch[i];
2789  Assert(cell->is_artificial() == false,
2790  ExcMessage("This function can not be called with cells that are "
2791  "not either locally owned or ghost cells."));
2792  local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
2793  cell->get_dof_indices(local_dof_indices);
2794  dofs_on_patch.insert(local_dof_indices.begin(),
2795  local_dof_indices.end());
2796  }
2797 
2798  Assert((dofs_on_patch.size() == count_dofs_on_patch<dim, spacedim>(patch)),
2799  ExcInternalError());
2800 
2801  // return a vector with the content of the set above. copying
2802  // also ensures that we retain sortedness as promised in the
2803  // documentation and as necessary to retain the block structure
2804  // also on the local system
2805  return std::vector<types::global_dof_index>(dofs_on_patch.begin(),
2806  dofs_on_patch.end());
2807  }
2808 
2809 
2810 } // end of namespace DoFTools
2811 
2812 
2813 
2814 // explicit instantiations
2815 
2816 #include "dof_tools.inst"
2817 
2818 
2819 
const std::vector< std::pair< size_type, number > > * get_constraint_entries(const size_type line_n) const
size_type n_constraints() const
unsigned int size() const
unsigned int n_selected_blocks(const unsigned int overall_number_of_blocks=numbers::invalid_unsigned_int) const
bool represents_n_components(const unsigned int n) const
bool represents_the_all_selected_mask() const
unsigned int size() const
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
unsigned int first_selected_component(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
std::vector< types::fe_index > get_active_fe_indices() const
types::global_dof_index n_boundary_dofs() const
const Triangulation< dim, spacedim > & get_triangulation() const
active_cell_iterator begin_active(const unsigned int level=0) const
bool has_hp_capabilities() const
types::global_dof_index n_dofs() const
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
cell_iterator begin(const unsigned int level=0) const
const IndexSet & locally_owned_dofs() const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
types::global_dof_index n_locally_owned_dofs() const
const std::vector< Point< spacedim > > & get_quadrature_points() const
const FEValues< dim, spacedim > & get_present_fe_values() const
unsigned int n_dofs_per_vertex() const
unsigned int n_dofs_per_cell() const
unsigned int n_dofs_per_line() const
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int n_blocks() const
unsigned int n_components() const
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
const ComponentMask & get_nonzero_components(const unsigned int i) const
bool is_primitive() const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
unsigned int component_to_block_index(const unsigned int component) const
unsigned int element_multiplicity(const unsigned int index) const
std::pair< unsigned int, unsigned int > face_system_to_component_index(const unsigned int index, const unsigned int face_no=0) const
unsigned int n_base_elements() const
size_type index_within_set(const size_type global_index) const
Definition: index_set.h:1863
size_type n_elements() const
Definition: index_set.h:1796
bool is_element(const size_type index) const
Definition: index_set.h:1756
void add_index(const size_type index)
Definition: index_set.h:1667
void fill_binary_vector(VectorType &vector) const
void subtract_set(const IndexSet &other)
Definition: index_set.cc:268
size_type nth_index_in_set(const size_type local_index) const
Definition: index_set.h:1844
void compress() const
Definition: index_set.h:1656
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition: index_set.h:1714
void reinit(const size_type m, const size_type n, const ArrayView< const unsigned int > &row_lengths)
void add(const size_type i, const size_type j)
virtual MPI_Comm get_communicator() const
virtual types::subdomain_id locally_owned_subdomain() const
unsigned int n_active_cells() const
unsigned int n_levels() const
unsigned int n_vertices() const
Definition: vector.h:109
size_type size() const
unsigned int size() const
Definition: collection.h:264
unsigned int max_dofs_per_face() const
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
unsigned int n_components() const
unsigned int max_dofs_per_cell() const
void push_back(const Quadrature< dim_in > &new_quadrature)
Definition: q_collection.h:216
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
unsigned int level
Definition: grid_out.cc:4618
unsigned int cell_index
Definition: grid_tools.cc:1191
IteratorRange< FilteredIterator< BaseIterator > > filter_iterators(IteratorRange< BaseIterator > i, const Predicate &p)
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define Assert(cond, exc)
Definition: exceptions.h:1586
static ::ExceptionBase & ExcNotImplemented()
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1759
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1886
#define AssertIndexRange(index, range)
Definition: exceptions.h:1827
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
#define AssertThrow(cond, exc)
Definition: exceptions.h:1675
@ update_quadrature_points
Transformed quadrature points.
void get_block_association(const DoFHandler< dim, spacedim > &dof, std::vector< unsigned char > &dofs_by_block)
Definition: dof_tools.cc:247
std::vector< unsigned char > get_local_component_association(const FiniteElement< dim, spacedim > &fe, const ComponentMask &component_mask)
Definition: dof_tools.cc:123
void resolve_components(const FiniteElement< dim, spacedim > &fe, const std::vector< unsigned char > &dofs_by_component, const std::vector< unsigned int > &target_component, const bool only_once, std::vector< types::global_dof_index > &dofs_per_component, unsigned int &component)
Definition: dof_tools.cc:1759
void get_component_association(const DoFHandler< dim, spacedim > &dof, const ComponentMask &component_mask, std::vector< unsigned char > &dofs_by_component)
Definition: dof_tools.cc:194
void get_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, std::vector< types::subdomain_id > &subdomain)
Definition: dof_tools.cc:1547
void extract_boundary_dofs(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< bool > &selected_dofs, const std::set< types::boundary_id > &boundary_ids={})
Definition: dof_tools.cc:545
IndexSet extract_dofs_with_support_contained_within(const DoFHandler< dim, spacedim > &dof_handler, const std::function< bool(const typename DoFHandler< dim, spacedim >::active_cell_iterator &)> &predicate, const AffineConstraints< number > &constraints=AffineConstraints< number >())
Definition: dof_tools.cc:796
IndexSet dof_indices_with_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain)
Definition: dof_tools.cc:1654
void write_gnuplot_dof_support_point_info(std::ostream &out, const std::map< types::global_dof_index, Point< spacedim >> &support_points)
Definition: dof_tools.cc:2334
std::vector< IndexSet > locally_owned_dofs_per_subdomain(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1376
IndexSet extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1136
void get_active_fe_indices(const DoFHandler< dim, spacedim > &dof_handler, std::vector< unsigned int > &active_fe_indices)
Definition: dof_tools.cc:1363
void make_cell_patches(SparsityPattern &block_list, const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, const std::vector< bool > &selected_dofs={}, const types::global_dof_index offset=0)
Definition: dof_tools.cc:2426
IndexSet extract_dofs(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask)
Definition: dof_tools.cc:386
IndexSet extract_locally_active_dofs(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1044
void convert_couplings_to_blocks(const DoFHandler< dim, spacedim > &dof_handler, const Table< 2, Coupling > &table_by_component, std::vector< Table< 2, Coupling >> &tables_by_block)
Definition: dof_tools.cc:2373
IndexSet extract_locally_relevant_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level)
Definition: dof_tools.cc:1186
void map_dofs_to_support_points(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, std::vector< Point< spacedim >> &support_points, const ComponentMask &mask=ComponentMask())
Definition: dof_tools.cc:2244
std::vector< IndexSet > locally_owned_dofs_per_component(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &components=ComponentMask())
Definition: dof_tools.cc:439
void extract_subdomain_dofs(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain_id, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:1013
IndexSet extract_locally_active_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level)
Definition: dof_tools.cc:1087
std::vector< unsigned int > make_vertex_patches(SparsityPattern &block_list, const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, const bool interior_dofs_only, const bool boundary_patches=false, const bool level_boundary_patches=false, const bool single_cell_patches=false, const bool invert_vertex_mapping=false)
Definition: dof_tools.cc:2574
std::vector< types::global_dof_index > count_dofs_per_fe_block(const DoFHandler< dim, spacedim > &dof, const std::vector< unsigned int > &target_block=std::vector< unsigned int >())
Definition: dof_tools.cc:1973
void extract_level_dofs(const unsigned int level, const DoFHandler< dim, spacedim > &dof, const ComponentMask &component_mask, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:473
void extract_dofs_with_support_on_boundary(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< bool > &selected_dofs, const std::set< types::boundary_id > &boundary_ids=std::set< types::boundary_id >())
Definition: dof_tools.cc:714
void distribute_cell_to_dof_vector(const DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &cell_data, Vector< double > &dof_data, const unsigned int component=0)
Definition: dof_tools.cc:302
std::vector< types::global_dof_index > count_dofs_per_fe_component(const DoFHandler< dim, spacedim > &dof_handler, const bool vector_valued_once=false, const std::vector< unsigned int > &target_component={})
Definition: dof_tools.cc:1881
void make_child_patches(SparsityPattern &block_list, const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, const bool interior_dofs_only, const bool boundary_dofs=false)
Definition: dof_tools.cc:2511
void map_dof_to_boundary_indices(const DoFHandler< dim, spacedim > &dof_handler, std::vector< types::global_dof_index > &mapping)
Definition: dof_tools.cc:2062
std::vector< IndexSet > locally_relevant_dofs_per_subdomain(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1471
void make_single_patch(SparsityPattern &block_list, const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, const bool interior_dofs_only=false)
Definition: dof_tools.cc:2469
void map_dofs_to_support_points(const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, std::map< types::global_dof_index, Point< spacedim >> &support_points, const ComponentMask &mask)
Definition: dof_tools.cc:2316
std::vector< types::global_dof_index > get_dofs_on_patch(const std::vector< typename DoFHandler< dim, spacedim >::active_cell_iterator > &patch)
Definition: dof_tools.cc:2775
unsigned int count_dofs_with_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain)
Definition: dof_tools.cc:1637
unsigned int count_dofs_on_patch(const std::vector< typename DoFHandler< dim, spacedim >::active_cell_iterator > &patch)
Definition: dof_tools.cc:2744
IndexSet extract_hanging_node_dofs(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1004
void extract_constant_modes(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< std::vector< bool >> &constant_modes)
Definition: dof_tools.cc:1247
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
std::vector< typename MeshType::active_cell_iterator > compute_active_cell_halo_layer(const MeshType &mesh, const std::function< bool(const typename MeshType::active_cell_iterator &)> &predicate)
std::enable_if_t< IsBlockVector< VectorType >::value, unsigned int > n_blocks(const VectorType &vector)
Definition: operators.h:50
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
VectorType::value_type * begin(VectorType &V)
VectorType::value_type * end(VectorType &V)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:150
std::string compress(const std::string &input)
Definition: utilities.cc:390
const types::boundary_id internal_face_boundary_id
Definition: types.h:277
const types::subdomain_id artificial_subdomain_id
Definition: types.h:315
const types::subdomain_id invalid_subdomain_id
Definition: types.h:298
static const unsigned int invalid_unsigned_int
Definition: types.h:213
const types::global_dof_index invalid_dof_index
Definition: types.h:233
unsigned int global_dof_index
Definition: types.h:82
unsigned int subdomain_id
Definition: types.h:44
unsigned short int fe_index
Definition: types.h:60
bool operator()(const Point< dim, Number > &lhs, const Point< dim, Number > &rhs) const
Definition: dof_tools.cc:81
const ::Triangulation< dim, spacedim > & tria
#define DEAL_II_DOF_INDEX_MPI_TYPE
Definition: types.h:99