Reference documentation for deal.II version GIT 4e92747666 2022-12-02 08:40: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  const ::DoFHandler<1, spacedim> &dof_handler)
906  {
907  // there are no hanging nodes in 1d
908  return IndexSet(dof_handler.n_dofs());
909  }
910 
911 
912  template <int spacedim>
913  IndexSet
915  const ::DoFHandler<2, spacedim> &dof_handler)
916  {
917  const unsigned int dim = 2;
918 
919  IndexSet selected_dofs(dof_handler.n_dofs());
920 
921  const FiniteElement<dim, spacedim> &fe = dof_handler.get_fe();
922 
923  // this function is similar to the make_sparsity_pattern function,
924  // see there for more information
925  for (const auto &cell : dof_handler.active_cell_iterators())
926  if (!cell->is_artificial())
927  {
928  for (const unsigned int face : cell->face_indices())
929  if (cell->face(face)->has_children())
930  {
931  const typename ::DoFHandler<dim,
932  spacedim>::line_iterator
933  line = cell->face(face);
934 
935  for (unsigned int dof = 0; dof != fe.n_dofs_per_vertex();
936  ++dof)
937  selected_dofs.add_index(
938  line->child(0)->vertex_dof_index(1, dof));
939 
940  for (unsigned int child = 0; child < 2; ++child)
941  {
942  if (cell->neighbor_child_on_subface(face, child)
943  ->is_artificial())
944  continue;
945  for (unsigned int dof = 0; dof != fe.n_dofs_per_line();
946  ++dof)
947  selected_dofs.add_index(
948  line->child(child)->dof_index(dof));
949  }
950  }
951  }
952 
953  selected_dofs.compress();
954  return selected_dofs;
955  }
956 
957 
958  template <int spacedim>
959  IndexSet
961  const ::DoFHandler<3, spacedim> &dof_handler)
962  {
963  const unsigned int dim = 3;
964 
965  IndexSet selected_dofs(dof_handler.n_dofs());
966  IndexSet unconstrained_dofs(dof_handler.n_dofs());
967 
968  const FiniteElement<dim, spacedim> &fe = dof_handler.get_fe();
969 
970  for (const auto &cell : dof_handler.active_cell_iterators())
971  if (!cell->is_artificial())
972  for (auto f : cell->face_indices())
973  {
974  const typename ::DoFHandler<dim, spacedim>::face_iterator
975  face = cell->face(f);
976  if (cell->face(f)->has_children())
977  {
978  for (unsigned int child = 0; child < 4; ++child)
979  if (!cell->neighbor_child_on_subface(f, child)
980  ->is_artificial())
981  {
982  // simply take all DoFs that live on this subface
983  std::vector<types::global_dof_index> ldi(
984  fe.n_dofs_per_face(f, child));
985  face->child(child)->get_dof_indices(ldi);
986  selected_dofs.add_indices(ldi.begin(), ldi.end());
987  }
988 
989  // and subtract (in the end) all the indices which a shared
990  // between this face and its subfaces
991  for (unsigned int vertex = 0; vertex < 4; ++vertex)
992  for (unsigned int dof = 0; dof != fe.n_dofs_per_vertex();
993  ++dof)
994  unconstrained_dofs.add_index(
995  face->vertex_dof_index(vertex, dof));
996  }
997  }
998  selected_dofs.subtract_set(unconstrained_dofs);
999  return selected_dofs;
1000  }
1001  } // namespace
1002  } // namespace internal
1003 
1004 
1005 
1006  template <int dim, int spacedim>
1007  IndexSet
1009  {
1010  return internal::extract_hanging_node_dofs(dof_handler);
1011  }
1012 
1013 
1014 
1015  template <int dim, int spacedim>
1016  void
1019  std::vector<bool> & selected_dofs)
1020  {
1021  Assert(selected_dofs.size() == dof_handler.n_dofs(),
1022  ExcDimensionMismatch(selected_dofs.size(), dof_handler.n_dofs()));
1023 
1024  // preset all values by false
1025  std::fill_n(selected_dofs.begin(), dof_handler.n_dofs(), false);
1026 
1027  std::vector<types::global_dof_index> local_dof_indices;
1028  local_dof_indices.reserve(
1029  dof_handler.get_fe_collection().max_dofs_per_cell());
1030 
1031  // this function is similar to the make_sparsity_pattern function, see
1032  // there for more information
1033  for (const auto &cell : dof_handler.active_cell_iterators())
1034  if (cell->subdomain_id() == subdomain_id)
1035  {
1036  const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
1037  local_dof_indices.resize(dofs_per_cell);
1038  cell->get_dof_indices(local_dof_indices);
1039  for (unsigned int i = 0; i < dofs_per_cell; ++i)
1040  selected_dofs[local_dof_indices[i]] = true;
1041  }
1042  }
1043 
1044 
1045 
1046  template <int dim, int spacedim>
1047  IndexSet
1049  {
1050  // collect all the locally owned dofs
1051  IndexSet dof_set = dof_handler.locally_owned_dofs();
1052 
1053  // add the DoF on the adjacent ghost cells to the IndexSet, cache them
1054  // in a set. need to check each dof manually because we can't be sure
1055  // that the dof range of locally_owned_dofs is really contiguous.
1056  std::vector<types::global_dof_index> dof_indices;
1057  std::set<types::global_dof_index> global_dof_indices;
1058 
1059  for (const auto &cell : dof_handler.active_cell_iterators() |
1061  {
1062  dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1063  cell->get_dof_indices(dof_indices);
1064 
1065  for (const types::global_dof_index dof_index : dof_indices)
1066  if (!dof_set.is_element(dof_index))
1067  global_dof_indices.insert(dof_index);
1068  }
1069 
1070  dof_set.add_indices(global_dof_indices.begin(), global_dof_indices.end());
1071 
1072  dof_set.compress();
1073 
1074  return dof_set;
1075  }
1076 
1077 
1078 
1079  template <int dim, int spacedim>
1080  void
1082  IndexSet & dof_set)
1083  {
1084  dof_set = extract_locally_active_dofs(dof_handler);
1085  }
1086 
1087 
1088 
1089  template <int dim, int spacedim>
1090  IndexSet
1092  const DoFHandler<dim, spacedim> &dof_handler,
1093  const unsigned int level)
1094  {
1095  // collect all the locally owned dofs
1096  IndexSet dof_set = dof_handler.locally_owned_mg_dofs(level);
1097 
1098  // add the DoF on the adjacent ghost cells to the IndexSet, cache them
1099  // in a set. need to check each dof manually because we can't be sure
1100  // that the dof range of locally_owned_dofs is really contiguous.
1101  std::vector<types::global_dof_index> dof_indices;
1102  std::set<types::global_dof_index> global_dof_indices;
1103 
1104  const auto filtered_iterators_range =
1107  for (const auto &cell : filtered_iterators_range)
1108  {
1109  dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1110  cell->get_mg_dof_indices(dof_indices);
1111 
1112  for (const types::global_dof_index dof_index : dof_indices)
1113  if (!dof_set.is_element(dof_index))
1114  global_dof_indices.insert(dof_index);
1115  }
1116 
1117  dof_set.add_indices(global_dof_indices.begin(), global_dof_indices.end());
1118 
1119  dof_set.compress();
1120 
1121  return dof_set;
1122  }
1123 
1124 
1125 
1126  template <int dim, int spacedim>
1127  void
1129  const DoFHandler<dim, spacedim> &dof_handler,
1130  IndexSet & dof_set,
1131  const unsigned int level)
1132  {
1133  dof_set = extract_locally_active_level_dofs(dof_handler, level);
1134  }
1135 
1136 
1137 
1138  template <int dim, int spacedim>
1139  IndexSet
1141  {
1142  // collect all the locally owned dofs
1143  IndexSet dof_set = dof_handler.locally_owned_dofs();
1144 
1145  // now add the DoF on the adjacent ghost cells to the IndexSet
1146 
1147  // Note: For certain meshes (in particular in 3D and with many
1148  // processors), it is really necessary to cache intermediate data. After
1149  // trying several objects such as std::set, a vector that is always kept
1150  // sorted, and a vector that is initially unsorted and sorted once at the
1151  // end, the latter has been identified to provide the best performance.
1152  // Martin Kronbichler
1153  std::vector<types::global_dof_index> dof_indices;
1154  std::vector<types::global_dof_index> dofs_on_ghosts;
1155 
1156  for (const auto &cell : dof_handler.active_cell_iterators())
1157  if (cell->is_ghost())
1158  {
1159  dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1160  cell->get_dof_indices(dof_indices);
1161  for (const auto dof_index : dof_indices)
1162  if (!dof_set.is_element(dof_index))
1163  dofs_on_ghosts.push_back(dof_index);
1164  }
1165 
1166  // sort, compress out duplicates, fill into index set
1167  std::sort(dofs_on_ghosts.begin(), dofs_on_ghosts.end());
1168  dof_set.add_indices(dofs_on_ghosts.begin(),
1169  std::unique(dofs_on_ghosts.begin(),
1170  dofs_on_ghosts.end()));
1171  dof_set.compress();
1172 
1173  return dof_set;
1174  }
1175 
1176 
1177 
1178  template <int dim, int spacedim>
1179  void
1181  IndexSet & dof_set)
1182  {
1183  dof_set = extract_locally_relevant_dofs(dof_handler);
1184  }
1185 
1186 
1187 
1188  template <int dim, int spacedim>
1189  IndexSet
1191  const DoFHandler<dim, spacedim> &dof_handler,
1192  const unsigned int level)
1193  {
1194  // collect all the locally owned dofs
1195  IndexSet dof_set = dof_handler.locally_owned_mg_dofs(level);
1196 
1197  // add the DoF on the adjacent ghost cells to the IndexSet
1198 
1199  // Note: For certain meshes (in particular in 3D and with many
1200  // processors), it is really necessary to cache intermediate data. After
1201  // trying several objects such as std::set, a vector that is always kept
1202  // sorted, and a vector that is initially unsorted and sorted once at the
1203  // end, the latter has been identified to provide the best performance.
1204  // Martin Kronbichler
1205  std::vector<types::global_dof_index> dof_indices;
1206  std::vector<types::global_dof_index> dofs_on_ghosts;
1207 
1208  for (const auto &cell : dof_handler.cell_iterators_on_level(level))
1209  {
1210  const types::subdomain_id id = cell->level_subdomain_id();
1211 
1212  // skip artificial and own cells (only look at ghost cells)
1213  if (id == dof_handler.get_triangulation().locally_owned_subdomain() ||
1215  continue;
1216 
1217  dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1218  cell->get_mg_dof_indices(dof_indices);
1219  for (const auto dof_index : dof_indices)
1220  if (!dof_set.is_element(dof_index))
1221  dofs_on_ghosts.push_back(dof_index);
1222  }
1223 
1224  // sort, compress out duplicates, fill into index set
1225  std::sort(dofs_on_ghosts.begin(), dofs_on_ghosts.end());
1226  dof_set.add_indices(dofs_on_ghosts.begin(),
1227  std::unique(dofs_on_ghosts.begin(),
1228  dofs_on_ghosts.end()));
1229 
1230  dof_set.compress();
1231 
1232  return dof_set;
1233  }
1234 
1235 
1236 
1237  template <int dim, int spacedim>
1238  void
1240  const DoFHandler<dim, spacedim> &dof_handler,
1241  const unsigned int level,
1242  IndexSet & dof_set)
1243  {
1244  dof_set = extract_locally_relevant_level_dofs(dof_handler, level);
1245  }
1246 
1247 
1248 
1249  template <int dim, int spacedim>
1250  void
1252  const ComponentMask & component_mask,
1253  std::vector<std::vector<bool>> & constant_modes)
1254  {
1255  // If there are no locally owned DoFs, return with an empty
1256  // constant_modes object:
1257  if (dof_handler.n_locally_owned_dofs() == 0)
1258  {
1259  constant_modes = std::vector<std::vector<bool>>(0);
1260  return;
1261  }
1262 
1263  const unsigned int n_components = dof_handler.get_fe(0).n_components();
1264  Assert(component_mask.represents_n_components(n_components),
1265  ExcDimensionMismatch(n_components, component_mask.size()));
1266 
1267  std::vector<unsigned char> dofs_by_component(
1268  dof_handler.n_locally_owned_dofs());
1270  component_mask,
1271  dofs_by_component);
1272  unsigned int n_selected_dofs = 0;
1273  for (unsigned int i = 0; i < n_components; ++i)
1274  if (component_mask[i] == true)
1275  n_selected_dofs +=
1276  std::count(dofs_by_component.begin(), dofs_by_component.end(), i);
1277 
1278  // Find local numbering within the selected components
1279  const IndexSet &locally_owned_dofs = dof_handler.locally_owned_dofs();
1280  std::vector<unsigned int> component_numbering(
1281  locally_owned_dofs.n_elements(), numbers::invalid_unsigned_int);
1282  for (unsigned int i = 0, count = 0; i < locally_owned_dofs.n_elements();
1283  ++i)
1284  if (component_mask[dofs_by_component[i]])
1285  component_numbering[i] = count++;
1286 
1287  // get the element constant modes and find a translation table between
1288  // index in the constant modes and the components.
1289  //
1290  // TODO: We might be able to extend this also for elements which do not
1291  // have the same constant modes, but that is messy...
1292  const ::hp::FECollection<dim, spacedim> &fe_collection =
1293  dof_handler.get_fe_collection();
1294  std::vector<Table<2, bool>> element_constant_modes;
1295  std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
1296  constant_mode_to_component_translation(n_components);
1297  {
1298  unsigned int n_constant_modes = 0;
1299  int first_non_empty_constant_mode = -1;
1300  for (unsigned int f = 0; f < fe_collection.size(); ++f)
1301  {
1302  std::pair<Table<2, bool>, std::vector<unsigned int>> data =
1303  fe_collection[f].get_constant_modes();
1304 
1305  // Store the index of the current element if it is the first that has
1306  // non-empty constant modes.
1307  if (first_non_empty_constant_mode < 0 && data.first.n_rows() > 0)
1308  {
1309  first_non_empty_constant_mode = f;
1310  // This is the first non-empty constant mode, so we figure out the
1311  // translation between index in the constant modes and the
1312  // components
1313  for (unsigned int i = 0; i < data.second.size(); ++i)
1314  if (component_mask[data.second[i]])
1315  constant_mode_to_component_translation[data.second[i]]
1316  .emplace_back(n_constant_modes++, i);
1317  }
1318 
1319  // Add the constant modes of this element to the list and assert that
1320  // there are as many constant modes as for the other elements (or zero
1321  // constant modes).
1322  element_constant_modes.push_back(data.first);
1323  Assert(
1324  element_constant_modes.back().n_rows() == 0 ||
1325  element_constant_modes.back().n_rows() ==
1326  element_constant_modes[first_non_empty_constant_mode].n_rows(),
1327  ExcInternalError());
1328  }
1329  AssertIndexRange(first_non_empty_constant_mode, fe_collection.size());
1330 
1331  // Now we know the number of constant modes and resize the return vector
1332  // accordingly
1333  constant_modes.clear();
1334  constant_modes.resize(n_constant_modes,
1335  std::vector<bool>(n_selected_dofs, false));
1336  }
1337 
1338  // Loop over all owned cells and ask the element for the constant modes
1339  std::vector<types::global_dof_index> dof_indices;
1340  for (const auto &cell : dof_handler.active_cell_iterators() |
1342  {
1343  dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1344  cell->get_dof_indices(dof_indices);
1345 
1346  for (unsigned int i = 0; i < dof_indices.size(); ++i)
1347  if (locally_owned_dofs.is_element(dof_indices[i]))
1348  {
1349  const unsigned int loc_index =
1350  locally_owned_dofs.index_within_set(dof_indices[i]);
1351  const unsigned int comp = dofs_by_component[loc_index];
1352  if (component_mask[comp])
1353  for (auto &indices :
1354  constant_mode_to_component_translation[comp])
1355  constant_modes[indices
1356  .first][component_numbering[loc_index]] =
1357  element_constant_modes[cell->active_fe_index()](
1358  indices.second, i);
1359  }
1360  }
1361  }
1362 
1363 
1364 
1365  template <int dim, int spacedim>
1366  void
1368  std::vector<unsigned int> & active_fe_indices)
1369  {
1370  AssertDimension(active_fe_indices.size(),
1371  dof_handler.get_triangulation().n_active_cells());
1372 
1373  std::vector<types::fe_index> indices = dof_handler.get_active_fe_indices();
1374 
1375  active_fe_indices.assign(indices.begin(), indices.end());
1376  }
1377 
1378  template <int dim, int spacedim>
1379  std::vector<IndexSet>
1381  {
1382  Assert(dof_handler.n_dofs() > 0,
1383  ExcMessage("The given DoFHandler has no DoFs."));
1384 
1385  // If the Triangulation is distributed, the only thing we can usefully
1386  // ask is for its locally owned subdomain
1387  Assert((dynamic_cast<
1389  &dof_handler.get_triangulation()) == nullptr),
1390  ExcMessage(
1391  "For parallel::distributed::Triangulation objects and "
1392  "associated DoF handler objects, asking for any information "
1393  "related to a subdomain other than the locally owned one does "
1394  "not make sense."));
1395 
1396  // The following is a random process (flip of a coin), thus should be called
1397  // once only.
1398  std::vector<::types::subdomain_id> subdomain_association(
1399  dof_handler.n_dofs());
1401  subdomain_association);
1402 
1403  // Figure out how many subdomain ids there are.
1404  //
1405  // if this is a parallel triangulation, then we can just ask the
1406  // triangulation for this. if this is a sequential triangulation, we loop
1407  // over all cells and take the largest subdomain_id value we find; the
1408  // number of subdomains is then the largest found value plus one. (we here
1409  // assume that all subdomain ids up to the largest are actually used; this
1410  // may not be true for a sequential triangulation where these values have
1411  // been set by hand and not in accordance with some MPI communicator; but
1412  // the function returns an array indexed starting at zero, so we need to
1413  // collect information for each subdomain index anyway, not just for the
1414  // used one.)
1415  const unsigned int n_subdomains =
1416  (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1417  &dof_handler.get_triangulation()) == nullptr ?
1418  [&dof_handler]() {
1419  unsigned int max_subdomain_id = 0;
1420  for (const auto &cell : dof_handler.active_cell_iterators())
1421  max_subdomain_id =
1422  std::max(max_subdomain_id, cell->subdomain_id());
1423  return max_subdomain_id + 1;
1424  }() :
1426  dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1427  &dof_handler.get_triangulation())
1428  ->get_communicator()));
1429  Assert(n_subdomains > *std::max_element(subdomain_association.begin(),
1430  subdomain_association.end()),
1431  ExcInternalError());
1432 
1433  std::vector<::IndexSet> index_sets(
1434  n_subdomains, ::IndexSet(dof_handler.n_dofs()));
1435 
1436  // loop over subdomain_association and populate IndexSet when a
1437  // change in subdomain ID is found
1438  ::types::global_dof_index i_min = 0;
1439  ::types::global_dof_index this_subdomain = subdomain_association[0];
1440 
1441  for (::types::global_dof_index index = 1;
1442  index < subdomain_association.size();
1443  ++index)
1444  {
1445  // found index different from the current one
1446  if (subdomain_association[index] != this_subdomain)
1447  {
1448  index_sets[this_subdomain].add_range(i_min, index);
1449  i_min = index;
1450  this_subdomain = subdomain_association[index];
1451  }
1452  }
1453 
1454  // the very last element is of different index
1455  if (i_min == subdomain_association.size() - 1)
1456  {
1457  index_sets[this_subdomain].add_index(i_min);
1458  }
1459 
1460  // otherwise there are at least two different indices
1461  else
1462  {
1463  index_sets[this_subdomain].add_range(i_min,
1464  subdomain_association.size());
1465  }
1466 
1467  for (unsigned int i = 0; i < n_subdomains; ++i)
1468  index_sets[i].compress();
1469 
1470  return index_sets;
1471  }
1472 
1473  template <int dim, int spacedim>
1474  std::vector<IndexSet>
1476  const DoFHandler<dim, spacedim> &dof_handler)
1477  {
1478  // If the Triangulation is distributed, the only thing we can usefully
1479  // ask is for its locally owned subdomain
1480  Assert((dynamic_cast<
1482  &dof_handler.get_triangulation()) == nullptr),
1483  ExcMessage(
1484  "For parallel::distributed::Triangulation objects and "
1485  "associated DoF handler objects, asking for any information "
1486  "related to a subdomain other than the locally owned one does "
1487  "not make sense."));
1488 
1489  // Collect all the locally owned DoFs
1490  // Note: Even though the distribution of DoFs by the
1491  // locally_owned_dofs_per_subdomain function is pseudo-random, we will
1492  // collect all the DoFs on the subdomain and its layer cell. Therefore, the
1493  // random nature of this function does not play a role in the extraction of
1494  // the locally relevant DoFs
1495  std::vector<IndexSet> dof_set =
1496  locally_owned_dofs_per_subdomain(dof_handler);
1497  const ::types::subdomain_id n_subdomains = dof_set.size();
1498 
1499  // Add the DoFs on the adjacent (equivalent ghost) cells to the IndexSet,
1500  // cache them in a set. Need to check each DoF manually because we can't
1501  // be sure that the DoF range of locally_owned_dofs is really contiguous.
1502  for (::types::subdomain_id subdomain_id = 0;
1503  subdomain_id < n_subdomains;
1504  ++subdomain_id)
1505  {
1506  // Extract the layer of cells around this subdomain
1507  std::function<bool(
1510  const std::vector<
1512  active_halo_layer =
1513  GridTools::compute_active_cell_halo_layer(dof_handler, predicate);
1514 
1515  // Extract DoFs associated with halo layer
1516  std::vector<types::global_dof_index> local_dof_indices;
1517  std::set<types::global_dof_index> subdomain_halo_global_dof_indices;
1518  for (typename std::vector<
1520  const_iterator it_cell = active_halo_layer.begin();
1521  it_cell != active_halo_layer.end();
1522  ++it_cell)
1523  {
1525  &cell = *it_cell;
1526  Assert(
1527  cell->subdomain_id() != subdomain_id,
1528  ExcMessage(
1529  "The subdomain ID of the halo cell should not match that of the vector entry."));
1530 
1531  local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1532  cell->get_dof_indices(local_dof_indices);
1533 
1534  for (const types::global_dof_index local_dof_index :
1535  local_dof_indices)
1536  subdomain_halo_global_dof_indices.insert(local_dof_index);
1537  }
1538 
1539  dof_set[subdomain_id].add_indices(
1540  subdomain_halo_global_dof_indices.begin(),
1541  subdomain_halo_global_dof_indices.end());
1542 
1543  dof_set[subdomain_id].compress();
1544  }
1545 
1546  return dof_set;
1547  }
1548 
1549  template <int dim, int spacedim>
1550  void
1552  const DoFHandler<dim, spacedim> & dof_handler,
1553  std::vector<types::subdomain_id> &subdomain_association)
1554  {
1555  // if the Triangulation is distributed, the only thing we can usefully
1556  // ask is for its locally owned subdomain
1557  Assert((dynamic_cast<
1559  &dof_handler.get_triangulation()) == nullptr),
1560  ExcMessage(
1561  "For parallel::distributed::Triangulation objects and "
1562  "associated DoF handler objects, asking for any subdomain other "
1563  "than the locally owned one does not make sense."));
1564 
1565  Assert(subdomain_association.size() == dof_handler.n_dofs(),
1566  ExcDimensionMismatch(subdomain_association.size(),
1567  dof_handler.n_dofs()));
1568 
1569  // catch an error that happened in some versions of the shared tria
1570  // distribute_dofs() function where we were trying to call this
1571  // function at a point in time when not all internal DoFHandler
1572  // structures were quite set up yet.
1573  Assert(dof_handler.n_dofs() > 0, ExcInternalError());
1574 
1575  // In case this function is executed with parallel::shared::Triangulation
1576  // with possibly artificial cells, we need to take "true" subdomain IDs
1577  // (i.e. without artificial cells). Otherwise we are good to use
1578  // subdomain_id as stored in cell->subdomain_id().
1579  std::vector<types::subdomain_id> cell_owners(
1580  dof_handler.get_triangulation().n_active_cells());
1582  (dynamic_cast<const parallel::shared::Triangulation<dim, spacedim> *>(
1583  &dof_handler.get_triangulation())))
1584  {
1585  cell_owners = tr->get_true_subdomain_ids_of_cells();
1586  Assert(tr->get_true_subdomain_ids_of_cells().size() ==
1587  tr->n_active_cells(),
1588  ExcInternalError());
1589  }
1590  else
1591  {
1592  for (const auto &cell : dof_handler.active_cell_iterators() |
1594  cell_owners[cell->active_cell_index()] = cell->subdomain_id();
1595  }
1596 
1597  // preset all values by an invalid value
1598  std::fill_n(subdomain_association.begin(),
1599  dof_handler.n_dofs(),
1601 
1602  std::vector<types::global_dof_index> local_dof_indices;
1603  local_dof_indices.reserve(
1604  dof_handler.get_fe_collection().max_dofs_per_cell());
1605 
1606  // loop over all cells and record which subdomain a DoF belongs to.
1607  // give to the smaller subdomain_id in case it is on an interface
1608  for (const auto &cell : dof_handler.active_cell_iterators())
1609  {
1611  cell_owners[cell->active_cell_index()];
1612  const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
1613  local_dof_indices.resize(dofs_per_cell);
1614  cell->get_dof_indices(local_dof_indices);
1615 
1616  // set subdomain ids. if dofs already have their values set then
1617  // they must be on partition interfaces. in that case assign them
1618  // to either the previous association or the current processor
1619  // with the smaller subdomain id.
1620  for (unsigned int i = 0; i < dofs_per_cell; ++i)
1621  if (subdomain_association[local_dof_indices[i]] ==
1623  subdomain_association[local_dof_indices[i]] = subdomain_id;
1624  else if (subdomain_association[local_dof_indices[i]] > subdomain_id)
1625  {
1626  subdomain_association[local_dof_indices[i]] = subdomain_id;
1627  }
1628  }
1629 
1630  Assert(std::find(subdomain_association.begin(),
1631  subdomain_association.end(),
1633  subdomain_association.end(),
1634  ExcInternalError());
1635  }
1636 
1637 
1638 
1639  template <int dim, int spacedim>
1640  unsigned int
1642  const DoFHandler<dim, spacedim> &dof_handler,
1643  const types::subdomain_id subdomain)
1644  {
1645  std::vector<types::subdomain_id> subdomain_association(
1646  dof_handler.n_dofs());
1647  get_subdomain_association(dof_handler, subdomain_association);
1648 
1649  return std::count(subdomain_association.begin(),
1650  subdomain_association.end(),
1651  subdomain);
1652  }
1653 
1654 
1655 
1656  template <int dim, int spacedim>
1657  IndexSet
1659  const DoFHandler<dim, spacedim> &dof_handler,
1660  const types::subdomain_id subdomain)
1661  {
1662  // If we have a distributed::Triangulation only allow locally_owned
1663  // subdomain.
1666  (subdomain ==
1667  dof_handler.get_triangulation().locally_owned_subdomain()),
1668  ExcMessage(
1669  "For parallel::distributed::Triangulation objects and "
1670  "associated DoF handler objects, asking for any subdomain other "
1671  "than the locally owned one does not make sense."));
1672 
1673  IndexSet index_set(dof_handler.n_dofs());
1674 
1675  std::vector<types::global_dof_index> local_dof_indices;
1676  local_dof_indices.reserve(
1677  dof_handler.get_fe_collection().max_dofs_per_cell());
1678 
1679  // first generate an unsorted list of all indices which we fill from
1680  // the back. could also insert them directly into the IndexSet, but
1681  // that inserts indices in the middle, which is an O(n^2) algorithm and
1682  // hence too expensive. Could also use std::set, but that is in general
1683  // more expensive than a vector
1684  std::vector<types::global_dof_index> subdomain_indices;
1685 
1686  for (const auto &cell : dof_handler.active_cell_iterators())
1687  if ((cell->is_artificial() == false) &&
1688  (cell->subdomain_id() == subdomain))
1689  {
1690  const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
1691  local_dof_indices.resize(dofs_per_cell);
1692  cell->get_dof_indices(local_dof_indices);
1693  subdomain_indices.insert(subdomain_indices.end(),
1694  local_dof_indices.begin(),
1695  local_dof_indices.end());
1696  }
1697  // sort indices and remove duplicates
1698  std::sort(subdomain_indices.begin(), subdomain_indices.end());
1699  subdomain_indices.erase(std::unique(subdomain_indices.begin(),
1700  subdomain_indices.end()),
1701  subdomain_indices.end());
1702 
1703  // insert into IndexSet
1704  index_set.add_indices(subdomain_indices.begin(), subdomain_indices.end());
1705  index_set.compress();
1706 
1707  return index_set;
1708  }
1709 
1710 
1711 
1712  template <int dim, int spacedim>
1713  void
1715  const DoFHandler<dim, spacedim> &dof_handler,
1716  const types::subdomain_id subdomain,
1717  std::vector<unsigned int> & n_dofs_on_subdomain)
1718  {
1719  Assert(n_dofs_on_subdomain.size() == dof_handler.get_fe(0).n_components(),
1720  ExcDimensionMismatch(n_dofs_on_subdomain.size(),
1721  dof_handler.get_fe(0).n_components()));
1722  std::fill(n_dofs_on_subdomain.begin(), n_dofs_on_subdomain.end(), 0);
1723 
1724  // Make sure there are at least some cells with this subdomain id
1725  Assert(std::any_of(
1726  dof_handler.begin_active(),
1728  dof_handler.end()},
1729  [subdomain](
1730  const typename DoFHandler<dim, spacedim>::cell_accessor &cell) {
1731  return cell.subdomain_id() == subdomain;
1732  }),
1733  ExcMessage("There are no cells for the given subdomain!"));
1734 
1735  std::vector<types::subdomain_id> subdomain_association(
1736  dof_handler.n_dofs());
1737  get_subdomain_association(dof_handler, subdomain_association);
1738 
1739  std::vector<unsigned char> component_association(dof_handler.n_dofs());
1741  std::vector<bool>(),
1742  component_association);
1743 
1744  for (unsigned int c = 0; c < dof_handler.get_fe(0).n_components(); ++c)
1745  {
1746  for (types::global_dof_index i = 0; i < dof_handler.n_dofs(); ++i)
1747  if ((subdomain_association[i] == subdomain) &&
1748  (component_association[i] == static_cast<unsigned char>(c)))
1749  ++n_dofs_on_subdomain[c];
1750  }
1751  }
1752 
1753 
1754 
1755  namespace internal
1756  {
1757  // TODO: why is this function so complicated? It would be nice to have
1758  // comments that explain why we can't just loop over all components and
1759  // count the entries in dofs_by_component that have this component's
1760  // index
1761  template <int dim, int spacedim>
1762  void
1764  const std::vector<unsigned char> & dofs_by_component,
1765  const std::vector<unsigned int> & target_component,
1766  const bool only_once,
1767  std::vector<types::global_dof_index> &dofs_per_component,
1768  unsigned int & component)
1769  {
1770  for (unsigned int b = 0; b < fe.n_base_elements(); ++b)
1771  {
1772  const FiniteElement<dim, spacedim> &base = fe.base_element(b);
1773  // Dimension of base element
1774  unsigned int d = base.n_components();
1775 
1776  for (unsigned int m = 0; m < fe.element_multiplicity(b); ++m)
1777  {
1778  if (base.n_base_elements() > 1)
1779  resolve_components(base,
1780  dofs_by_component,
1781  target_component,
1782  only_once,
1783  dofs_per_component,
1784  component);
1785  else
1786  {
1787  for (unsigned int dd = 0; dd < d; ++dd, ++component)
1788  dofs_per_component[target_component[component]] +=
1789  std::count(dofs_by_component.begin(),
1790  dofs_by_component.end(),
1791  component);
1792 
1793  // if we have non-primitive FEs and want all components
1794  // to show the number of dofs, need to copy the result to
1795  // those components
1796  if (!base.is_primitive() && !only_once)
1797  for (unsigned int dd = 1; dd < d; ++dd)
1798  dofs_per_component[target_component[component - d + dd]] =
1799  dofs_per_component[target_component[component - d]];
1800  }
1801  }
1802  }
1803  }
1804 
1805 
1806  template <int dim, int spacedim>
1807  void
1809  const std::vector<unsigned char> & dofs_by_component,
1810  const std::vector<unsigned int> & target_component,
1811  const bool only_once,
1812  std::vector<types::global_dof_index> &dofs_per_component,
1813  unsigned int & component)
1814  {
1815  // assert that all elements in the collection have the same structure
1816  // (base elements and multiplicity, components per base element) and
1817  // then simply call the function above
1818  for (unsigned int fe = 1; fe < fe_collection.size(); ++fe)
1819  {
1820  Assert(fe_collection[fe].n_components() ==
1821  fe_collection[0].n_components(),
1822  ExcNotImplemented());
1823  Assert(fe_collection[fe].n_base_elements() ==
1824  fe_collection[0].n_base_elements(),
1825  ExcNotImplemented());
1826  for (unsigned int b = 0; b < fe_collection[0].n_base_elements(); ++b)
1827  {
1828  Assert(fe_collection[fe].base_element(b).n_components() ==
1829  fe_collection[0].base_element(b).n_components(),
1830  ExcNotImplemented());
1831  Assert(fe_collection[fe].base_element(b).n_base_elements() ==
1832  fe_collection[0].base_element(b).n_base_elements(),
1833  ExcNotImplemented());
1834  }
1835  }
1836 
1837  resolve_components(fe_collection[0],
1838  dofs_by_component,
1839  target_component,
1840  only_once,
1841  dofs_per_component,
1842  component);
1843  }
1844  } // namespace internal
1845 
1846 
1847 
1848  namespace internal
1849  {
1850  namespace
1851  {
1855  template <int dim, int spacedim>
1856  bool
1857  all_elements_are_primitive(const FiniteElement<dim, spacedim> &fe)
1858  {
1859  return fe.is_primitive();
1860  }
1861 
1862 
1867  template <int dim, int spacedim>
1868  bool
1869  all_elements_are_primitive(
1870  const ::hp::FECollection<dim, spacedim> &fe_collection)
1871  {
1872  for (unsigned int i = 0; i < fe_collection.size(); ++i)
1873  if (fe_collection[i].is_primitive() == false)
1874  return false;
1875 
1876  return true;
1877  }
1878  } // namespace
1879  } // namespace internal
1880 
1881 
1882 
1883  template <int dim, int spacedim>
1884  std::vector<types::global_dof_index>
1886  const DoFHandler<dim, spacedim> &dof_handler,
1887  const bool only_once,
1888  const std::vector<unsigned int> &target_component_)
1889  {
1890  const unsigned int n_components = dof_handler.get_fe(0).n_components();
1891 
1892  // If the empty vector was given as default argument, set up this
1893  // vector as identity.
1894  std::vector<unsigned int> target_component = target_component_;
1895  if (target_component.size() == 0)
1896  {
1897  target_component.resize(n_components);
1898  for (unsigned int i = 0; i < n_components; ++i)
1899  target_component[i] = i;
1900  }
1901  else
1902  Assert(target_component.size() == n_components,
1903  ExcDimensionMismatch(target_component.size(), n_components));
1904 
1905 
1906  const unsigned int max_component =
1907  *std::max_element(target_component.begin(), target_component.end());
1908  const unsigned int n_target_components = max_component + 1;
1909 
1910  std::vector<types::global_dof_index> dofs_per_component(
1911  n_target_components, types::global_dof_index(0));
1912 
1913  // special case for only one component. treat this first since it does
1914  // not require any computations
1915  if (n_components == 1)
1916  {
1917  dofs_per_component[0] = dof_handler.n_locally_owned_dofs();
1918  return dofs_per_component;
1919  }
1920 
1921 
1922  // otherwise determine the number of dofs in each component separately.
1923  // do so in parallel
1924  std::vector<unsigned char> dofs_by_component(
1925  dof_handler.n_locally_owned_dofs());
1927  ComponentMask(),
1928  dofs_by_component);
1929 
1930  // next count what we got
1931  unsigned int component = 0;
1933  dofs_by_component,
1934  target_component,
1935  only_once,
1936  dofs_per_component,
1937  component);
1938  Assert(n_components == component, ExcInternalError());
1939 
1940  // finally sanity check. this is only valid if the finite element is
1941  // actually primitive, so exclude other elements from this
1942  Assert((internal::all_elements_are_primitive(
1943  dof_handler.get_fe_collection()) == false) ||
1944  (std::accumulate(dofs_per_component.begin(),
1945  dofs_per_component.end(),
1947  dof_handler.n_locally_owned_dofs()),
1948  ExcInternalError());
1949 
1950  // reduce information from all CPUs
1951 #ifdef DEAL_II_WITH_MPI
1952 
1954  (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1955  &dof_handler.get_triangulation())))
1956  {
1957  std::vector<types::global_dof_index> local_dof_count =
1958  dofs_per_component;
1959 
1960  const int ierr = MPI_Allreduce(local_dof_count.data(),
1961  dofs_per_component.data(),
1962  n_target_components,
1964  MPI_SUM,
1965  tria->get_communicator());
1966  AssertThrowMPI(ierr);
1967  }
1968 #endif
1969 
1970  return dofs_per_component;
1971  }
1972 
1973 
1974 
1975  template <int dim, int spacedim>
1976  std::vector<types::global_dof_index>
1978  const std::vector<unsigned int> &target_block_)
1979  {
1980  const ::hp::FECollection<dim, spacedim> &fe_collection =
1981  dof_handler.get_fe_collection();
1982  Assert(fe_collection.size() < 256, ExcNotImplemented());
1983 
1984  // If the empty vector for target_block(e.g., as default argument), then
1985  // set up this vector as identity. We do this set up with the first
1986  // element of the collection, but the whole thing can only work if
1987  // all elements have the same number of blocks anyway -- so check
1988  // that right after
1989  const unsigned int n_blocks = fe_collection[0].n_blocks();
1990 
1991  std::vector<unsigned int> target_block = target_block_;
1992  if (target_block.size() == 0)
1993  {
1994  target_block.resize(fe_collection[0].n_blocks());
1995  for (unsigned int i = 0; i < n_blocks; ++i)
1996  target_block[i] = i;
1997  }
1998  else
1999  Assert(target_block.size() == n_blocks,
2000  ExcDimensionMismatch(target_block.size(), n_blocks));
2001  for (unsigned int f = 1; f < fe_collection.size(); ++f)
2002  Assert(fe_collection[0].n_blocks() == fe_collection[f].n_blocks(),
2003  ExcMessage("This function can only work if all elements in a "
2004  "collection have the same number of blocks."));
2005 
2006  // special case for only one block. treat this first since it does
2007  // not require any computations
2008  if (n_blocks == 1)
2009  {
2010  std::vector<types::global_dof_index> dofs_per_block(1);
2011  dofs_per_block[0] = dof_handler.n_dofs();
2012  return dofs_per_block;
2013  }
2014 
2015  // Otherwise set up the right-sized object and start working
2016  const unsigned int max_block =
2017  *std::max_element(target_block.begin(), target_block.end());
2018  const unsigned int n_target_blocks = max_block + 1;
2019 
2020  std::vector<types::global_dof_index> dofs_per_block(n_target_blocks);
2021 
2022  // Loop over the elements of the collection, but really only consider
2023  // the last element (see #9271)
2024  for (unsigned int this_fe = fe_collection.size() - 1;
2025  this_fe < fe_collection.size();
2026  ++this_fe)
2027  {
2028  const FiniteElement<dim, spacedim> &fe = fe_collection[this_fe];
2029 
2030  std::vector<unsigned char> dofs_by_block(
2031  dof_handler.n_locally_owned_dofs());
2032  internal::get_block_association(dof_handler, dofs_by_block);
2033 
2034  // next count what we got
2035  for (unsigned int block = 0; block < fe.n_blocks(); ++block)
2036  dofs_per_block[target_block[block]] +=
2037  std::count(dofs_by_block.begin(), dofs_by_block.end(), block);
2038 
2039 #ifdef DEAL_II_WITH_MPI
2040  // if we are working on a parallel mesh, we now need to collect
2041  // this information from all processors
2043  (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
2044  &dof_handler.get_triangulation())))
2045  {
2046  std::vector<types::global_dof_index> local_dof_count =
2047  dofs_per_block;
2048  const int ierr = MPI_Allreduce(local_dof_count.data(),
2049  dofs_per_block.data(),
2050  n_target_blocks,
2052  MPI_SUM,
2053  tria->get_communicator());
2054  AssertThrowMPI(ierr);
2055  }
2056 #endif
2057  }
2058 
2059  return dofs_per_block;
2060  }
2061 
2062 
2063 
2064  template <int dim, int spacedim>
2065  void
2067  std::vector<types::global_dof_index> &mapping)
2068  {
2069  mapping.clear();
2070  mapping.insert(mapping.end(),
2071  dof_handler.n_dofs(),
2073 
2074  std::vector<types::global_dof_index> dofs_on_face;
2075  dofs_on_face.reserve(dof_handler.get_fe_collection().max_dofs_per_face());
2076  types::global_dof_index next_boundary_index = 0;
2077 
2078  // now loop over all cells and check whether their faces are at the
2079  // boundary. note that we need not take special care of single lines
2080  // being at the boundary (using @p{cell->has_boundary_lines}), since we
2081  // do not support boundaries of dimension dim-2, and so every isolated
2082  // boundary line is also part of a boundary face which we will be
2083  // visiting sooner or later
2084  for (const auto &cell : dof_handler.active_cell_iterators())
2085  for (const unsigned int f : cell->face_indices())
2086  if (cell->at_boundary(f))
2087  {
2088  const unsigned int dofs_per_face =
2089  cell->get_fe().n_dofs_per_face(f);
2090  dofs_on_face.resize(dofs_per_face);
2091  cell->face(f)->get_dof_indices(dofs_on_face,
2092  cell->active_fe_index());
2093  for (unsigned int i = 0; i < dofs_per_face; ++i)
2094  if (mapping[dofs_on_face[i]] == numbers::invalid_dof_index)
2095  mapping[dofs_on_face[i]] = next_boundary_index++;
2096  }
2097 
2098  AssertDimension(next_boundary_index, dof_handler.n_boundary_dofs());
2099  }
2100 
2101 
2102 
2103  template <int dim, int spacedim>
2104  void
2106  const std::set<types::boundary_id> &boundary_ids,
2107  std::vector<types::global_dof_index> &mapping)
2108  {
2109  Assert(boundary_ids.find(numbers::internal_face_boundary_id) ==
2110  boundary_ids.end(),
2112 
2113  mapping.clear();
2114  mapping.insert(mapping.end(),
2115  dof_handler.n_dofs(),
2117 
2118  // return if there is nothing to do
2119  if (boundary_ids.size() == 0)
2120  return;
2121 
2122  std::vector<types::global_dof_index> dofs_on_face;
2123  dofs_on_face.reserve(dof_handler.get_fe_collection().max_dofs_per_face());
2124  types::global_dof_index next_boundary_index = 0;
2125 
2126  for (const auto &cell : dof_handler.active_cell_iterators())
2127  for (const unsigned int f : cell->face_indices())
2128  if (boundary_ids.find(cell->face(f)->boundary_id()) !=
2129  boundary_ids.end())
2130  {
2131  const unsigned int dofs_per_face =
2132  cell->get_fe().n_dofs_per_face(f);
2133  dofs_on_face.resize(dofs_per_face);
2134  cell->face(f)->get_dof_indices(dofs_on_face,
2135  cell->active_fe_index());
2136  for (unsigned int i = 0; i < dofs_per_face; ++i)
2137  if (mapping[dofs_on_face[i]] == numbers::invalid_dof_index)
2138  mapping[dofs_on_face[i]] = next_boundary_index++;
2139  }
2140 
2141  AssertDimension(next_boundary_index,
2142  dof_handler.n_boundary_dofs(boundary_ids));
2143  }
2144 
2145  namespace internal
2146  {
2147  namespace
2148  {
2149  template <int dim, int spacedim>
2150  void
2152  const hp::MappingCollection<dim, spacedim> & mapping,
2153  const DoFHandler<dim, spacedim> & dof_handler,
2154  std::map<types::global_dof_index, Point<spacedim>> &support_points,
2155  const ComponentMask & in_mask)
2156  {
2157  const hp::FECollection<dim, spacedim> &fe_collection =
2158  dof_handler.get_fe_collection();
2159  hp::QCollection<dim> q_coll_dummy;
2160 
2161  for (unsigned int fe_index = 0; fe_index < fe_collection.size();
2162  ++fe_index)
2163  {
2164  // check whether every FE in the collection has support points
2165  Assert((fe_collection[fe_index].n_dofs_per_cell() == 0) ||
2166  (fe_collection[fe_index].has_support_points()),
2168  q_coll_dummy.push_back(Quadrature<dim>(
2169  fe_collection[fe_index].get_unit_support_points()));
2170  }
2171 
2172  // Take care of components
2173  const ComponentMask mask =
2174  (in_mask.size() == 0 ?
2175  ComponentMask(fe_collection.n_components(), true) :
2176  in_mask);
2177 
2178  // Now loop over all cells and enquire the support points on each
2179  // of these. we use dummy quadrature formulas where the quadrature
2180  // points are located at the unit support points to enquire the
2181  // location of the support points in real space.
2182  //
2183  // The weights of the quadrature rule have been set to invalid
2184  // values by the used constructor.
2185  hp::FEValues<dim, spacedim> hp_fe_values(mapping,
2186  fe_collection,
2187  q_coll_dummy,
2189 
2190  std::vector<types::global_dof_index> local_dof_indices;
2191  for (const auto &cell : dof_handler.active_cell_iterators())
2192  // only work on locally relevant cells
2193  if (cell->is_artificial() == false)
2194  {
2195  hp_fe_values.reinit(cell);
2196  const FEValues<dim, spacedim> &fe_values =
2197  hp_fe_values.get_present_fe_values();
2198 
2199  local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
2200  cell->get_dof_indices(local_dof_indices);
2201 
2202  const std::vector<Point<spacedim>> &points =
2203  fe_values.get_quadrature_points();
2204  for (unsigned int i = 0; i < cell->get_fe().n_dofs_per_cell();
2205  ++i)
2206  {
2207  const unsigned int dof_comp =
2208  cell->get_fe().system_to_component_index(i).first;
2209 
2210  // insert the values into the map if it is a valid component
2211  if (mask[dof_comp])
2212  support_points[local_dof_indices[i]] = points[i];
2213  }
2214  }
2215  }
2216 
2217 
2218  template <int dim, int spacedim>
2219  void
2221  const hp::MappingCollection<dim, spacedim> &mapping,
2222  const DoFHandler<dim, spacedim> & dof_handler,
2223  std::vector<Point<spacedim>> & support_points,
2224  const ComponentMask & mask)
2225  {
2226  // get the data in the form of the map as above
2227  std::map<types::global_dof_index, Point<spacedim>> x_support_points;
2229  dof_handler,
2230  x_support_points,
2231  mask);
2232 
2233  // now convert from the map to the linear vector. make sure every
2234  // entry really appeared in the map
2235  for (types::global_dof_index i = 0; i < dof_handler.n_dofs(); ++i)
2236  {
2237  Assert(x_support_points.find(i) != x_support_points.end(),
2238  ExcInternalError());
2239 
2240  support_points[i] = x_support_points[i];
2241  }
2242  }
2243  } // namespace
2244  } // namespace internal
2245 
2246  template <int dim, int spacedim>
2247  void
2249  const DoFHandler<dim, spacedim> &dof_handler,
2250  std::vector<Point<spacedim>> & support_points,
2251  const ComponentMask & mask)
2252  {
2253  AssertDimension(support_points.size(), dof_handler.n_dofs());
2254  Assert((dynamic_cast<
2256  &dof_handler.get_triangulation()) == nullptr),
2257  ExcMessage(
2258  "This function can not be used with distributed triangulations. "
2259  "See the documentation for more information."));
2260 
2261  // Let the internal function do all the work, just make sure that it
2262  // gets a MappingCollection
2263  const hp::MappingCollection<dim, spacedim> mapping_collection(mapping);
2264 
2265  internal::map_dofs_to_support_points(mapping_collection,
2266  dof_handler,
2267  support_points,
2268  mask);
2269  }
2270 
2271 
2272  template <int dim, int spacedim>
2273  void
2275  const hp::MappingCollection<dim, spacedim> &mapping,
2276  const DoFHandler<dim, spacedim> & dof_handler,
2277  std::vector<Point<spacedim>> & support_points,
2278  const ComponentMask & mask)
2279  {
2280  AssertDimension(support_points.size(), dof_handler.n_dofs());
2281  Assert((dynamic_cast<
2283  &dof_handler.get_triangulation()) == nullptr),
2284  ExcMessage(
2285  "This function can not be used with distributed triangulations. "
2286  "See the documentation for more information."));
2287 
2288  // Let the internal function do all the work, just make sure that it
2289  // gets a MappingCollection
2291  dof_handler,
2292  support_points,
2293  mask);
2294  }
2295 
2296 
2297  template <int dim, int spacedim>
2298  void
2300  const Mapping<dim, spacedim> & mapping,
2301  const DoFHandler<dim, spacedim> & dof_handler,
2302  std::map<types::global_dof_index, Point<spacedim>> &support_points,
2303  const ComponentMask & mask)
2304  {
2305  support_points.clear();
2306 
2307  // Let the internal function do all the work, just make sure that it
2308  // gets a MappingCollection
2309  const hp::MappingCollection<dim, spacedim> mapping_collection(mapping);
2310 
2311  internal::map_dofs_to_support_points(mapping_collection,
2312  dof_handler,
2313  support_points,
2314  mask);
2315  }
2316 
2317 
2318  template <int dim, int spacedim>
2319  void
2321  const hp::MappingCollection<dim, spacedim> & mapping,
2322  const DoFHandler<dim, spacedim> & dof_handler,
2323  std::map<types::global_dof_index, Point<spacedim>> &support_points,
2324  const ComponentMask & mask)
2325  {
2326  support_points.clear();
2327 
2328  // Let the internal function do all the work, just make sure that it
2329  // gets a MappingCollection
2331  dof_handler,
2332  support_points,
2333  mask);
2334  }
2335 
2336  template <int spacedim>
2337  void
2339  std::ostream & out,
2340  const std::map<types::global_dof_index, Point<spacedim>> &support_points)
2341  {
2342  AssertThrow(out.fail() == false, ExcIO());
2343 
2344  std::map<Point<spacedim>,
2345  std::vector<types::global_dof_index>,
2347  point_map;
2348 
2349  // convert to map point -> list of DoFs
2350  for (const auto &it : support_points)
2351  {
2352  std::vector<types::global_dof_index> &v = point_map[it.second];
2353  v.push_back(it.first);
2354  }
2355 
2356  // print the newly created map:
2357  for (const auto &it : point_map)
2358  {
2359  out << it.first << " \"";
2360  const std::vector<types::global_dof_index> &v = it.second;
2361  for (unsigned int i = 0; i < v.size(); ++i)
2362  {
2363  if (i > 0)
2364  out << ", ";
2365  out << v[i];
2366  }
2367 
2368  out << "\"\n";
2369  }
2370 
2371  out << std::flush;
2372  }
2373 
2374 
2375  template <int dim, int spacedim>
2376  void
2378  const Table<2, Coupling> & table,
2379  std::vector<Table<2, Coupling>> &tables_by_block)
2380  {
2381  if (dof_handler.has_hp_capabilities() == false)
2382  {
2383  const FiniteElement<dim, spacedim> &fe = dof_handler.get_fe();
2384  const unsigned int nb = fe.n_blocks();
2385 
2386  tables_by_block.resize(1);
2387  tables_by_block[0].reinit(nb, nb);
2388  tables_by_block[0].fill(none);
2389 
2390  for (unsigned int i = 0; i < fe.n_components(); ++i)
2391  {
2392  const unsigned int ib = fe.component_to_block_index(i);
2393  for (unsigned int j = 0; j < fe.n_components(); ++j)
2394  {
2395  const unsigned int jb = fe.component_to_block_index(j);
2396  tables_by_block[0](ib, jb) |= table(i, j);
2397  }
2398  }
2399  }
2400  else
2401  {
2402  const hp::FECollection<dim> &fe_collection =
2403  dof_handler.get_fe_collection();
2404  tables_by_block.resize(fe_collection.size());
2405 
2406  for (unsigned int f = 0; f < fe_collection.size(); ++f)
2407  {
2408  const FiniteElement<dim, spacedim> &fe = fe_collection[f];
2409 
2410  const unsigned int nb = fe.n_blocks();
2411  tables_by_block[f].reinit(nb, nb);
2412  tables_by_block[f].fill(none);
2413  for (unsigned int i = 0; i < fe.n_components(); ++i)
2414  {
2415  const unsigned int ib = fe.component_to_block_index(i);
2416  for (unsigned int j = 0; j < fe.n_components(); ++j)
2417  {
2418  const unsigned int jb = fe.component_to_block_index(j);
2419  tables_by_block[f](ib, jb) |= table(i, j);
2420  }
2421  }
2422  }
2423  }
2424  }
2425 
2426 
2427 
2428  template <int dim, int spacedim>
2429  void
2431  const DoFHandler<dim, spacedim> &dof_handler,
2432  const unsigned int level,
2433  const std::vector<bool> & selected_dofs,
2434  const types::global_dof_index offset)
2435  {
2436  std::vector<types::global_dof_index> indices;
2437 
2438  unsigned int i = 0;
2439 
2440  for (const auto &cell : dof_handler.cell_iterators_on_level(level))
2441  if (cell->is_locally_owned_on_level())
2442  ++i;
2443  block_list.reinit(i,
2444  dof_handler.n_dofs(),
2445  dof_handler.get_fe().n_dofs_per_cell());
2446  i = 0;
2447  for (const auto &cell : dof_handler.cell_iterators_on_level(level))
2448  if (cell->is_locally_owned_on_level())
2449  {
2450  indices.resize(cell->get_fe().n_dofs_per_cell());
2451  cell->get_mg_dof_indices(indices);
2452 
2453  if (selected_dofs.size() != 0)
2454  AssertDimension(indices.size(), selected_dofs.size());
2455 
2456  for (types::global_dof_index j = 0; j < indices.size(); ++j)
2457  {
2458  if (selected_dofs.size() == 0)
2459  block_list.add(i, indices[j] - offset);
2460  else
2461  {
2462  if (selected_dofs[j])
2463  block_list.add(i, indices[j] - offset);
2464  }
2465  }
2466  ++i;
2467  }
2468  }
2469 
2470 
2471  template <int dim, int spacedim>
2472  void
2474  const DoFHandler<dim, spacedim> &dof_handler,
2475  const unsigned int level,
2476  const bool interior_only)
2477  {
2478  const FiniteElement<dim> &fe = dof_handler.get_fe();
2479  block_list.reinit(1, dof_handler.n_dofs(level), dof_handler.n_dofs(level));
2480 
2481  std::vector<types::global_dof_index> indices;
2482  std::vector<bool> exclude;
2483 
2484  for (const auto &cell : dof_handler.cell_iterators_on_level(level))
2485  {
2486  indices.resize(cell->get_fe().n_dofs_per_cell());
2487  cell->get_mg_dof_indices(indices);
2488 
2489  if (interior_only)
2490  {
2491  // Exclude degrees of freedom on faces opposite to the vertex
2492  exclude.resize(fe.n_dofs_per_cell());
2493  std::fill(exclude.begin(), exclude.end(), false);
2494 
2495  for (const unsigned int face : cell->face_indices())
2496  if (cell->at_boundary(face) ||
2497  cell->neighbor(face)->level() != cell->level())
2498  for (unsigned int i = 0; i < fe.n_dofs_per_face(face); ++i)
2499  exclude[fe.face_to_cell_index(i, face)] = true;
2500  for (types::global_dof_index j = 0; j < indices.size(); ++j)
2501  if (!exclude[j])
2502  block_list.add(0, indices[j]);
2503  }
2504  else
2505  {
2506  for (const auto index : indices)
2507  block_list.add(0, index);
2508  }
2509  }
2510  }
2511 
2512 
2513  template <int dim, int spacedim>
2514  void
2516  const DoFHandler<dim, spacedim> &dof_handler,
2517  const unsigned int level,
2518  const bool interior_dofs_only,
2519  const bool boundary_dofs)
2520  {
2521  Assert(level > 0 && level < dof_handler.get_triangulation().n_levels(),
2522  ExcIndexRange(level, 1, dof_handler.get_triangulation().n_levels()));
2523 
2524  std::vector<types::global_dof_index> indices;
2525  std::vector<bool> exclude;
2526 
2527  unsigned int block = 0;
2528  for (const auto &pcell : dof_handler.cell_iterators_on_level(level - 1))
2529  {
2530  if (pcell->is_active())
2531  continue;
2532 
2533  for (unsigned int child = 0; child < pcell->n_children(); ++child)
2534  {
2535  const auto cell = pcell->child(child);
2536 
2537  // For hp, only this line here would have to be replaced.
2538  const FiniteElement<dim> &fe = dof_handler.get_fe();
2539  const unsigned int n_dofs = fe.n_dofs_per_cell();
2540  indices.resize(n_dofs);
2541  exclude.resize(n_dofs);
2542  std::fill(exclude.begin(), exclude.end(), false);
2543  cell->get_mg_dof_indices(indices);
2544 
2545  if (interior_dofs_only)
2546  {
2547  // Eliminate dofs on faces of the child which are on faces
2548  // of the parent
2549  for (unsigned int d = 0; d < dim; ++d)
2550  {
2551  const unsigned int face =
2553  for (unsigned int i = 0; i < fe.n_dofs_per_face(face); ++i)
2554  exclude[fe.face_to_cell_index(i, face)] = true;
2555  }
2556 
2557  // Now remove all degrees of freedom on the domain boundary
2558  // from the exclusion list
2559  if (boundary_dofs)
2560  for (const unsigned int face :
2562  if (cell->at_boundary(face))
2563  for (unsigned int i = 0; i < fe.n_dofs_per_face(face);
2564  ++i)
2565  exclude[fe.face_to_cell_index(i, face)] = false;
2566  }
2567 
2568  for (unsigned int i = 0; i < n_dofs; ++i)
2569  if (!exclude[i])
2570  block_list.add(block, indices[i]);
2571  }
2572  ++block;
2573  }
2574  }
2575 
2576  template <int dim, int spacedim>
2577  std::vector<unsigned int>
2579  const DoFHandler<dim, spacedim> &dof_handler,
2580  const unsigned int level,
2581  const bool interior_only,
2582  const bool boundary_patches,
2583  const bool level_boundary_patches,
2584  const bool single_cell_patches,
2585  const bool invert_vertex_mapping)
2586  {
2587  const unsigned int n_blocks = dof_handler.get_fe().n_blocks();
2588  BlockMask exclude_boundary_dofs = BlockMask(n_blocks, interior_only);
2589  return make_vertex_patches(block_list,
2590  dof_handler,
2591  level,
2592  exclude_boundary_dofs,
2593  boundary_patches,
2594  level_boundary_patches,
2595  single_cell_patches,
2596  invert_vertex_mapping);
2597  }
2598 
2599  template <int dim, int spacedim>
2600  std::vector<unsigned int>
2602  const DoFHandler<dim, spacedim> &dof_handler,
2603  const unsigned int level,
2604  const BlockMask & exclude_boundary_dofs,
2605  const bool boundary_patches,
2606  const bool level_boundary_patches,
2607  const bool single_cell_patches,
2608  const bool invert_vertex_mapping)
2609  {
2610  // Vector mapping from vertex index in the triangulation to consecutive
2611  // block indices on this level The number of cells at a vertex
2612  std::vector<unsigned int> vertex_cell_count(
2613  dof_handler.get_triangulation().n_vertices(), 0);
2614 
2615  // Is a vertex at the boundary?
2616  std::vector<bool> vertex_boundary(
2617  dof_handler.get_triangulation().n_vertices(), false);
2618 
2619  std::vector<unsigned int> vertex_mapping(
2620  dof_handler.get_triangulation().n_vertices(),
2622 
2623  // Estimate for the number of dofs at this point
2624  std::vector<unsigned int> vertex_dof_count(
2625  dof_handler.get_triangulation().n_vertices(), 0);
2626 
2627  // Identify all vertices active on this level and remember some data
2628  // about them
2629  for (const auto &cell : dof_handler.cell_iterators_on_level(level))
2630  for (const unsigned int v : cell->vertex_indices())
2631  {
2632  const unsigned int vg = cell->vertex_index(v);
2633  vertex_dof_count[vg] += cell->get_fe().n_dofs_per_cell();
2634  ++vertex_cell_count[vg];
2635  for (unsigned int d = 0; d < dim; ++d)
2636  {
2637  const unsigned int face = GeometryInfo<dim>::vertex_to_face[v][d];
2638  if (cell->at_boundary(face))
2639  vertex_boundary[vg] = true;
2640  else if ((!level_boundary_patches) &&
2641  (cell->neighbor(face)->level() !=
2642  static_cast<int>(level)))
2643  vertex_boundary[vg] = true;
2644  }
2645  }
2646  // From now on, only vertices with positive dof count are "in".
2647 
2648  // Remove vertices at boundaries or in corners
2649  for (unsigned int vg = 0; vg < vertex_dof_count.size(); ++vg)
2650  if ((!single_cell_patches && vertex_cell_count[vg] < 2) ||
2651  (!boundary_patches && vertex_boundary[vg]))
2652  vertex_dof_count[vg] = 0;
2653 
2654  // Create a mapping from all vertices to the ones used here
2655  unsigned int n_vertex_count = 0;
2656  for (unsigned int vg = 0; vg < vertex_mapping.size(); ++vg)
2657  if (vertex_dof_count[vg] != 0)
2658  vertex_mapping[vg] = n_vertex_count++;
2659 
2660  // Compactify dof count
2661  for (unsigned int vg = 0; vg < vertex_mapping.size(); ++vg)
2662  if (vertex_dof_count[vg] != 0)
2663  vertex_dof_count[vertex_mapping[vg]] = vertex_dof_count[vg];
2664 
2665  // Now that we have all the data, we reduce it to the part we actually
2666  // want
2667  vertex_dof_count.resize(n_vertex_count);
2668 
2669  // At this point, the list of patches is ready. Now we enter the dofs
2670  // into the sparsity pattern.
2671  block_list.reinit(vertex_dof_count.size(),
2672  dof_handler.n_dofs(level),
2673  vertex_dof_count);
2674 
2675  std::vector<types::global_dof_index> indices;
2676  std::vector<bool> exclude;
2677 
2678  for (const auto &cell : dof_handler.cell_iterators_on_level(level))
2679  {
2680  const FiniteElement<dim> &fe = cell->get_fe();
2681  indices.resize(fe.n_dofs_per_cell());
2682  cell->get_mg_dof_indices(indices);
2683 
2684  for (const unsigned int v : cell->vertex_indices())
2685  {
2686  const unsigned int vg = cell->vertex_index(v);
2687  const unsigned int block = vertex_mapping[vg];
2688  if (block == numbers::invalid_unsigned_int)
2689  continue;
2690 
2691  // Collect excluded dofs for some block(s) if boundary dofs
2692  // for a block are decided to be excluded
2693  if (exclude_boundary_dofs.size() == 0 ||
2694  exclude_boundary_dofs.n_selected_blocks() != 0)
2695  {
2696  // Exclude degrees of freedom on faces opposite to the
2697  // vertex
2698  exclude.resize(fe.n_dofs_per_cell());
2699  std::fill(exclude.begin(), exclude.end(), false);
2700 
2701  for (unsigned int d = 0; d < dim; ++d)
2702  {
2703  const unsigned int a_face =
2705  const unsigned int face =
2707  for (unsigned int i = 0; i < fe.n_dofs_per_face(face); ++i)
2708  {
2709  // For each dof, get the block it is in and decide to
2710  // exclude it or not
2711  if (exclude_boundary_dofs[fe.system_to_block_index(
2712  fe.face_to_cell_index(
2713  i, face))
2714  .first] == true)
2715  exclude[fe.face_to_cell_index(i, face)] = true;
2716  }
2717  }
2718  for (unsigned int j = 0; j < indices.size(); ++j)
2719  if (!exclude[j])
2720  block_list.add(block, indices[j]);
2721  }
2722  else
2723  {
2724  for (const auto index : indices)
2725  block_list.add(block, index);
2726  }
2727  }
2728  }
2729 
2730  if (invert_vertex_mapping)
2731  {
2732  // Compress vertex mapping
2733  unsigned int n_vertex_count = 0;
2734  for (unsigned int vg = 0; vg < vertex_mapping.size(); ++vg)
2735  if (vertex_mapping[vg] != numbers::invalid_unsigned_int)
2736  vertex_mapping[n_vertex_count++] = vg;
2737 
2738  // Now we reduce it to the part we actually want
2739  vertex_mapping.resize(n_vertex_count);
2740  }
2741 
2742  return vertex_mapping;
2743  }
2744 
2745 
2746  template <int dim, int spacedim>
2747  unsigned int
2749  const std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
2750  &patch)
2751  {
2752  std::set<types::global_dof_index> dofs_on_patch;
2753  std::vector<types::global_dof_index> local_dof_indices;
2754 
2755  // loop over the cells in the patch and get the DoFs on each.
2756  // add all of them to a std::set which automatically makes sure
2757  // all duplicates are ignored
2758  for (unsigned int i = 0; i < patch.size(); ++i)
2759  {
2761  patch[i];
2762  Assert(cell->is_artificial() == false,
2763  ExcMessage("This function can not be called with cells that are "
2764  "not either locally owned or ghost cells."));
2765  local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
2766  cell->get_dof_indices(local_dof_indices);
2767  dofs_on_patch.insert(local_dof_indices.begin(),
2768  local_dof_indices.end());
2769  }
2770 
2771  // now return the number of DoFs (duplicates were ignored)
2772  return dofs_on_patch.size();
2773  }
2774 
2775 
2776 
2777  template <int dim, int spacedim>
2778  std::vector<types::global_dof_index>
2780  const std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
2781  &patch)
2782  {
2783  std::set<types::global_dof_index> dofs_on_patch;
2784  std::vector<types::global_dof_index> local_dof_indices;
2785 
2786  // loop over the cells in the patch and get the DoFs on each.
2787  // add all of them to a std::set which automatically makes sure
2788  // all duplicates are ignored
2789  for (unsigned int i = 0; i < patch.size(); ++i)
2790  {
2792  patch[i];
2793  Assert(cell->is_artificial() == false,
2794  ExcMessage("This function can not be called with cells that are "
2795  "not either locally owned or ghost cells."));
2796  local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
2797  cell->get_dof_indices(local_dof_indices);
2798  dofs_on_patch.insert(local_dof_indices.begin(),
2799  local_dof_indices.end());
2800  }
2801 
2802  Assert((dofs_on_patch.size() == count_dofs_on_patch<dim, spacedim>(patch)),
2803  ExcInternalError());
2804 
2805  // return a vector with the content of the set above. copying
2806  // also ensures that we retain sortedness as promised in the
2807  // documentation and as necessary to retain the block structure
2808  // also on the local system
2809  return std::vector<types::global_dof_index>(dofs_on_patch.begin(),
2810  dofs_on_patch.end());
2811  }
2812 
2813 
2814 } // end of namespace DoFTools
2815 
2816 
2817 
2818 // explicit instantiations
2819 
2820 #include "dof_tools.inst"
2821 
2822 
2823 
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:1841
size_type n_elements() const
Definition: index_set.h:1774
bool is_element(const size_type index) const
Definition: index_set.h:1734
void add_index(const size_type index)
Definition: index_set.h:1645
void fill_binary_vector(VectorType &vector) const
void subtract_set(const IndexSet &other)
Definition: index_set.cc:266
size_type nth_index_in_set(const size_type local_index) const
Definition: index_set.h:1822
void compress() const
Definition: index_set.h:1634
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition: index_set.h:1692
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:458
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:459
unsigned int level
Definition: grid_out.cc:4608
unsigned int cell_index
Definition: grid_tools.cc:1177
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:1501
static ::ExceptionBase & ExcNotImplemented()
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1695
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1818
#define AssertIndexRange(index, range)
Definition: exceptions.h:1760
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:1611
@ 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:1763
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:1551
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:1658
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:2338
std::vector< IndexSet > locally_owned_dofs_per_subdomain(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1380
IndexSet extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1140
void get_active_fe_indices(const DoFHandler< dim, spacedim > &dof_handler, std::vector< unsigned int > &active_fe_indices)
Definition: dof_tools.cc:1367
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:2430
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:1048
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:2377
IndexSet extract_locally_relevant_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level)
Definition: dof_tools.cc:1190
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:2248
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:1017
IndexSet extract_locally_active_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level)
Definition: dof_tools.cc:1091
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:2578
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:1977
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:1885
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:2515
void map_dof_to_boundary_indices(const DoFHandler< dim, spacedim > &dof_handler, std::vector< types::global_dof_index > &mapping)
Definition: dof_tools.cc:2066
std::vector< IndexSet > locally_relevant_dofs_per_subdomain(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1475
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:2473
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:2320
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:2779
unsigned int count_dofs_with_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain)
Definition: dof_tools.cc:1641
unsigned int count_dofs_on_patch(const std::vector< typename DoFHandler< dim, spacedim >::active_cell_iterator > &patch)
Definition: dof_tools.cc:2748
IndexSet extract_hanging_node_dofs(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1008
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:1251
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:146
std::string compress(const std::string &input)
Definition: utilities.cc:392
const types::boundary_id internal_face_boundary_id
Definition: types.h:270
const types::subdomain_id artificial_subdomain_id
Definition: types.h:308
const types::subdomain_id invalid_subdomain_id
Definition: types.h:291
static const unsigned int invalid_unsigned_int
Definition: types.h:206
const types::global_dof_index invalid_dof_index
Definition: types.h:226
unsigned int global_dof_index
Definition: types.h:81
unsigned int subdomain_id
Definition: types.h:43
unsigned short int fe_index
Definition: types.h:59
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:91