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