Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
dof_tools.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2019 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_dof_tools_h
17 #define dealii_dof_tools_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/base/exceptions.h>
23 #include <deal.II/base/index_set.h>
24 #include <deal.II/base/point.h>
25 
26 #include <deal.II/dofs/dof_handler.h>
27 
28 #include <deal.II/fe/component_mask.h>
29 
30 #include <deal.II/hp/dof_handler.h>
31 
32 #include <deal.II/lac/affine_constraints.h>
33 
34 #include <map>
35 #include <ostream>
36 #include <set>
37 #include <vector>
38 
39 
40 DEAL_II_NAMESPACE_OPEN
41 
42 class BlockMask;
43 template <int dim, typename RangeNumberType>
44 class Function;
45 template <int dim, int spacedim>
46 class FiniteElement;
47 namespace hp
48 {
49  template <int dim, int spacedim>
51  template <int dim, int spacedim>
52  class FECollection;
53 } // namespace hp
54 template <class MeshType>
56 template <int dim, int spacedim>
57 class Mapping;
58 class SparsityPattern;
59 template <int dim, class T>
60 class Table;
61 template <typename Number>
62 class Vector;
63 
64 namespace GridTools
65 {
66  template <typename CellIterator>
67  struct PeriodicFacePair;
68 }
69 
70 
190 namespace DoFTools
191 {
203  enum Coupling
204  {
221  };
222 
237  template <int dim, int spacedim>
238  void
240  const Table<2, Coupling> &table_by_component,
241  std::vector<Table<2, Coupling>> &tables_by_block);
242 
252  template <int dim, int spacedim>
253  void
255  const Table<2, Coupling> &table_by_component,
256  std::vector<Table<2, Coupling>> &tables_by_block);
257 
263  template <int dim, int spacedim>
267  const Table<2, Coupling> & component_couplings);
268 
276  template <int dim, int spacedim>
277  std::vector<Table<2, Coupling>>
280  const Table<2, Coupling> & component_couplings);
299  template <int dim, int spacedim>
300  DEAL_II_DEPRECATED unsigned int
302 
310  template <int dim, int spacedim>
311  DEAL_II_DEPRECATED unsigned int
313 
314 
325  template <int dim, int spacedim>
326  DEAL_II_DEPRECATED unsigned int
328 
339  template <int dim, int spacedim>
340  DEAL_II_DEPRECATED unsigned int
342 
353  template <int dim, int spacedim>
354  DEAL_II_DEPRECATED unsigned int
356 
367  template <int dim, int spacedim>
368  DEAL_II_DEPRECATED unsigned int
370 
382  template <int dim, int spacedim>
383  DEAL_II_DEPRECATED unsigned int
385 
397  template <int dim, int spacedim>
398  DEAL_II_DEPRECATED unsigned int
400 
412  template <int dim, int spacedim>
413  DEAL_II_DEPRECATED bool
415 
427  template <int dim, int spacedim>
428  DEAL_II_DEPRECATED bool
430 
556  template <typename DoFHandlerType,
557  typename SparsityPatternType,
558  typename number = double>
559  void
561  const DoFHandlerType & dof_handler,
562  SparsityPatternType & sparsity_pattern,
564  const bool keep_constrained_dofs = true,
566 
632  template <typename DoFHandlerType,
633  typename SparsityPatternType,
634  typename number = double>
635  void
637  const DoFHandlerType & dof_handler,
638  const Table<2, Coupling> & coupling,
639  SparsityPatternType & sparsity_pattern,
641  const bool keep_constrained_dofs = true,
643 
664  template <typename DoFHandlerType, typename SparsityPatternType>
665  void
666  make_sparsity_pattern(const DoFHandlerType &dof_row,
667  const DoFHandlerType &dof_col,
668  SparsityPatternType & sparsity);
669 
715  template <typename DoFHandlerType, typename SparsityPatternType>
716  void
717  make_flux_sparsity_pattern(const DoFHandlerType &dof_handler,
718  SparsityPatternType & sparsity_pattern);
719 
728  template <typename DoFHandlerType,
729  typename SparsityPatternType,
730  typename number>
731  void
733  const DoFHandlerType & dof_handler,
734  SparsityPatternType & sparsity_pattern,
735  const AffineConstraints<number> &constraints,
736  const bool keep_constrained_dofs = true,
738 
739 
759  template <typename DoFHandlerType, typename SparsityPatternType>
760  void
762  const DoFHandlerType & dof,
763  SparsityPatternType & sparsity,
764  const Table<2, Coupling> &cell_integrals_mask,
765  const Table<2, Coupling> &face_integrals_mask,
767 
768 
777  template <typename DoFHandlerType,
778  typename SparsityPatternType,
779  typename number>
780  void
781  make_flux_sparsity_pattern(const DoFHandlerType & dof,
782  SparsityPatternType & sparsity,
783  const AffineConstraints<number> &constraints,
784  const bool keep_constrained_dofs,
785  const Table<2, Coupling> &couplings,
786  const Table<2, Coupling> &face_couplings,
787  const types::subdomain_id subdomain_id);
788 
798  template <typename DoFHandlerType, typename SparsityPatternType>
799  void
801  const DoFHandlerType & dof,
802  const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
803  SparsityPatternType & sparsity_pattern);
804 
822  template <typename DoFHandlerType,
823  typename SparsityPatternType,
824  typename number>
825  void
827  const DoFHandlerType &dof,
828  const std::map<types::boundary_id,
830  & boundary_ids,
831  const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
832  SparsityPatternType & sparsity);
833 
884  template <typename DoFHandlerType, typename number>
885  void
886  make_hanging_node_constraints(const DoFHandlerType & dof_handler,
887  AffineConstraints<number> &constraints);
888 
956  template <int dim, int spacedim>
957  void
959  const DoFHandler<dim, spacedim> & coarse_grid,
960  const unsigned int coarse_component,
961  const DoFHandler<dim, spacedim> & fine_grid,
962  const unsigned int fine_component,
963  const InterGridMap<DoFHandler<dim, spacedim>> &coarse_to_fine_grid_map,
964  AffineConstraints<double> & constraints);
965 
966 
983  template <int dim, int spacedim>
984  void
986  const DoFHandler<dim, spacedim> & coarse_grid,
987  const unsigned int coarse_component,
988  const DoFHandler<dim, spacedim> & fine_grid,
989  const unsigned int fine_component,
990  const InterGridMap<DoFHandler<dim, spacedim>> &coarse_to_fine_grid_map,
991  std::vector<std::map<types::global_dof_index, float>>
992  &transfer_representation);
993 
1156  template <typename FaceIterator>
1157  void
1159  const FaceIterator & face_1,
1160  const typename identity<FaceIterator>::type &face_2,
1161  AffineConstraints<double> & constraints,
1162  const ComponentMask & component_mask = ComponentMask(),
1163  const bool face_orientation = true,
1164  const bool face_flip = false,
1165  const bool face_rotation = false,
1166  const FullMatrix<double> & matrix = FullMatrix<double>(),
1167  const std::vector<unsigned int> &first_vector_components =
1168  std::vector<unsigned int>());
1169 
1170 
1171 
1193  template <typename DoFHandlerType>
1194  void
1196  const std::vector<
1198  & periodic_faces,
1199  AffineConstraints<double> & constraints,
1200  const ComponentMask & component_mask = ComponentMask(),
1201  const std::vector<unsigned int> &first_vector_components =
1202  std::vector<unsigned int>());
1203 
1204 
1205 
1236  template <typename DoFHandlerType>
1237  void
1239  const DoFHandlerType & dof_handler,
1240  const types::boundary_id b_id1,
1241  const types::boundary_id b_id2,
1242  const int direction,
1243  AffineConstraints<double> &constraints,
1244  const ComponentMask & component_mask = ComponentMask());
1245 
1246 
1247 
1272  template <typename DoFHandlerType>
1273  void
1275  const DoFHandlerType & dof_handler,
1276  const types::boundary_id b_id,
1277  const int direction,
1278  AffineConstraints<double> &constraints,
1279  const ComponentMask & component_mask = ComponentMask());
1280 
1305  template <int dim, int spacedim>
1306  void
1308  std::vector<bool> & selected_dofs);
1309 
1314  template <int dim, int spacedim>
1315  IndexSet
1317 
1351  template <int dim, int spacedim>
1352  void
1353  extract_dofs(const DoFHandler<dim, spacedim> &dof_handler,
1354  const ComponentMask & component_mask,
1355  std::vector<bool> & selected_dofs);
1356 
1360  template <int dim, int spacedim>
1361  void
1362  extract_dofs(const hp::DoFHandler<dim, spacedim> &dof_handler,
1363  const ComponentMask & component_mask,
1364  std::vector<bool> & selected_dofs);
1365 
1390  template <int dim, int spacedim>
1391  void
1392  extract_dofs(const DoFHandler<dim, spacedim> &dof_handler,
1393  const BlockMask & block_mask,
1394  std::vector<bool> & selected_dofs);
1395 
1399  template <int dim, int spacedim>
1400  void
1401  extract_dofs(const hp::DoFHandler<dim, spacedim> &dof_handler,
1402  const BlockMask & block_mask,
1403  std::vector<bool> & selected_dofs);
1404 
1409  template <typename DoFHandlerType>
1410  void
1411  extract_level_dofs(const unsigned int level,
1412  const DoFHandlerType &dof,
1413  const ComponentMask & component_mask,
1414  std::vector<bool> & selected_dofs);
1415 
1420  template <typename DoFHandlerType>
1421  void
1422  extract_level_dofs(const unsigned int level,
1423  const DoFHandlerType &dof,
1424  const BlockMask & component_mask,
1425  std::vector<bool> & selected_dofs);
1426 
1479  template <typename DoFHandlerType>
1480  void
1481  extract_boundary_dofs(const DoFHandlerType & dof_handler,
1482  const ComponentMask & component_mask,
1483  std::vector<bool> & selected_dofs,
1484  const std::set<types::boundary_id> &boundary_ids =
1485  std::set<types::boundary_id>());
1486 
1517  template <typename DoFHandlerType>
1518  void
1519  extract_boundary_dofs(const DoFHandlerType & dof_handler,
1520  const ComponentMask & component_mask,
1521  IndexSet & selected_dofs,
1522  const std::set<types::boundary_id> &boundary_ids =
1523  std::set<types::boundary_id>());
1524 
1541  template <typename DoFHandlerType>
1542  void
1544  const DoFHandlerType & dof_handler,
1545  const ComponentMask & component_mask,
1546  std::vector<bool> & selected_dofs,
1547  const std::set<types::boundary_id> &boundary_ids =
1548  std::set<types::boundary_id>());
1549 
1577  template <typename DoFHandlerType, typename number = double>
1578  IndexSet
1580  const DoFHandlerType &dof_handler,
1581  const std::function<
1582  bool(const typename DoFHandlerType::active_cell_iterator &)> &predicate,
1583  const AffineConstraints<number> &constraints = AffineConstraints<number>());
1584 
1618  template <typename DoFHandlerType>
1619  void
1620  extract_constant_modes(const DoFHandlerType & dof_handler,
1621  const ComponentMask & component_mask,
1622  std::vector<std::vector<bool>> &constant_modes);
1624 
1638  template <typename DoFHandlerType>
1639  void
1640  extract_subdomain_dofs(const DoFHandlerType & dof_handler,
1641  const types::subdomain_id subdomain_id,
1642  std::vector<bool> & selected_dofs);
1643 
1644 
1654  template <typename DoFHandlerType>
1655  DEAL_II_DEPRECATED void
1656  extract_locally_owned_dofs(const DoFHandlerType &dof_handler,
1657  IndexSet & dof_set);
1658 
1659 
1674  template <typename DoFHandlerType>
1675  void
1676  extract_locally_active_dofs(const DoFHandlerType &dof_handler,
1677  IndexSet & dof_set);
1678 
1688  template <typename DoFHandlerType>
1689  void
1690  extract_locally_relevant_dofs(const DoFHandlerType &dof_handler,
1691  IndexSet & dof_set);
1692 
1709  template <typename DoFHandlerType>
1710  std::vector<IndexSet>
1711  locally_owned_dofs_per_subdomain(const DoFHandlerType &dof_handler);
1712 
1729  template <typename DoFHandlerType>
1730  std::vector<IndexSet>
1731  locally_relevant_dofs_per_subdomain(const DoFHandlerType &dof_handler);
1732 
1733 
1738  template <typename DoFHandlerType>
1739  void
1740  extract_locally_relevant_level_dofs(const DoFHandlerType &dof_handler,
1741  const unsigned int level,
1742  IndexSet & dof_set);
1743 
1744 
1775  template <typename DoFHandlerType>
1776  void
1777  get_subdomain_association(const DoFHandlerType & dof_handler,
1778  std::vector<types::subdomain_id> &subdomain);
1779 
1805  template <typename DoFHandlerType>
1806  unsigned int
1807  count_dofs_with_subdomain_association(const DoFHandlerType & dof_handler,
1808  const types::subdomain_id subdomain);
1809 
1829  template <typename DoFHandlerType>
1830  void
1832  const DoFHandlerType & dof_handler,
1833  const types::subdomain_id subdomain,
1834  std::vector<unsigned int> &n_dofs_on_subdomain);
1835 
1857  template <typename DoFHandlerType>
1858  IndexSet
1859  dof_indices_with_subdomain_association(const DoFHandlerType & dof_handler,
1860  const types::subdomain_id subdomain);
1861  // @}
1871 
1923  template <typename DoFHandlerType>
1924  std::vector<types::global_dof_index>
1926  const std::vector<typename DoFHandlerType::active_cell_iterator> &patch);
1927 
1948  template <int dim, int spacedim>
1949  void
1950  make_cell_patches(SparsityPattern & block_list,
1951  const DoFHandler<dim, spacedim> &dof_handler,
1952  const unsigned int level,
1953  const std::vector<bool> & selected_dofs = {},
1954  const types::global_dof_index offset = 0);
1955 
2007  template <typename DoFHandlerType>
2008  std::vector<unsigned int>
2009  make_vertex_patches(SparsityPattern & block_list,
2010  const DoFHandlerType &dof_handler,
2011  const unsigned int level,
2012  const bool interior_dofs_only,
2013  const bool boundary_patches = false,
2014  const bool level_boundary_patches = false,
2015  const bool single_cell_patches = false,
2016  const bool invert_vertex_mapping = false);
2017 
2032  template <typename DoFHandlerType>
2033  std::vector<unsigned int>
2034  make_vertex_patches(SparsityPattern & block_list,
2035  const DoFHandlerType &dof_handler,
2036  const unsigned int level,
2037  const BlockMask & exclude_boundary_dofs = BlockMask(),
2038  const bool boundary_patches = false,
2039  const bool level_boundary_patches = false,
2040  const bool single_cell_patches = false,
2041  const bool invert_vertex_mapping = false);
2042 
2080  template <typename DoFHandlerType>
2081  void
2082  make_child_patches(SparsityPattern & block_list,
2083  const DoFHandlerType &dof_handler,
2084  const unsigned int level,
2085  const bool interior_dofs_only,
2086  const bool boundary_dofs = false);
2087 
2107  template <typename DoFHandlerType>
2108  void
2109  make_single_patch(SparsityPattern & block_list,
2110  const DoFHandlerType &dof_handler,
2111  const unsigned int level,
2112  const bool interior_dofs_only = false);
2113 
2156  template <typename DoFHandlerType>
2157  void
2159  const DoFHandlerType & dof_handler,
2160  std::vector<types::global_dof_index> &dofs_per_component,
2161  const bool vector_valued_once = false,
2162  std::vector<unsigned int> target_component = {});
2163 
2180  template <typename DoFHandlerType>
2181  void
2182  count_dofs_per_block(const DoFHandlerType & dof,
2183  std::vector<types::global_dof_index> &dofs_per_block,
2184  const std::vector<unsigned int> & target_block =
2185  std::vector<unsigned int>());
2186 
2197  template <typename DoFHandlerType>
2198  void
2199  get_active_fe_indices(const DoFHandlerType & dof_handler,
2200  std::vector<unsigned int> &active_fe_indices);
2201 
2236  template <typename DoFHandlerType>
2237  unsigned int
2239  const std::vector<typename DoFHandlerType::active_cell_iterator> &patch);
2240 
2263  template <typename DoFHandlerType>
2264  void
2265  map_dof_to_boundary_indices(const DoFHandlerType & dof_handler,
2266  std::vector<types::global_dof_index> &mapping);
2267 
2278  template <typename DoFHandlerType>
2279  void
2280  map_dof_to_boundary_indices(const DoFHandlerType & dof_handler,
2281  const std::set<types::boundary_id> &boundary_ids,
2282  std::vector<types::global_dof_index> &mapping);
2283 
2301  template <int dim, int spacedim>
2302  void
2304  const DoFHandler<dim, spacedim> &dof_handler,
2305  std::vector<Point<spacedim>> & support_points);
2306 
2311  template <int dim, int spacedim>
2312  void
2314  const ::hp::MappingCollection<dim, spacedim> &mapping,
2315  const hp::DoFHandler<dim, spacedim> & dof_handler,
2316  std::vector<Point<spacedim>> & support_points);
2317 
2345  template <int dim, int spacedim>
2346  void
2348  const Mapping<dim, spacedim> & mapping,
2349  const DoFHandler<dim, spacedim> & dof_handler,
2350  std::map<types::global_dof_index, Point<spacedim>> &support_points);
2351 
2355  template <int dim, int spacedim>
2356  void
2358  const ::hp::MappingCollection<dim, spacedim> &mapping,
2359  const hp::DoFHandler<dim, spacedim> & dof_handler,
2360  std::map<types::global_dof_index, Point<spacedim>> &support_points);
2361 
2362 
2383  template <typename DoFHandlerType, class Comp>
2384  void
2387  & mapping,
2388  const DoFHandlerType &dof_handler,
2391  Comp> & point_to_index_map);
2429  template <typename DoFHandlerType, typename Number>
2430  void
2431  distribute_cell_to_dof_vector(const DoFHandlerType &dof_handler,
2432  const Vector<Number> &cell_data,
2433  Vector<double> & dof_data,
2434  const unsigned int component = 0);
2435 
2436 
2513  template <int spacedim>
2514  void
2516  std::ostream & out,
2517  const std::map<types::global_dof_index, Point<spacedim>> &support_points);
2518 
2519 
2560  template <int dim,
2561  int spacedim,
2562  template <int, int> class DoFHandlerType,
2563  typename number>
2564  void
2566  const DoFHandlerType<dim, spacedim> &dof,
2567  const types::boundary_id boundary_id,
2568  AffineConstraints<number> & zero_boundary_constraints,
2569  const ComponentMask & component_mask = ComponentMask());
2570 
2581  template <int dim,
2582  int spacedim,
2583  template <int, int> class DoFHandlerType,
2584  typename number>
2585  void
2587  const DoFHandlerType<dim, spacedim> &dof,
2588  AffineConstraints<number> & zero_boundary_constraints,
2589  const ComponentMask & component_mask = ComponentMask());
2590 
2635 } // namespace DoFTools
2636 
2637 
2638 
2639 /* ------------------------- inline functions -------------- */
2640 
2641 #ifndef DOXYGEN
2642 
2643 namespace DoFTools
2644 {
2650  inline Coupling
2651  operator|=(Coupling &c1, const Coupling c2)
2652  {
2653  if (c2 == always)
2654  c1 = always;
2655  else if (c1 != always && c2 == nonzero)
2656  return c1 = nonzero;
2657  return c1;
2658  }
2659 
2660 
2666  inline Coupling
2667  operator|(const Coupling c1, const Coupling c2)
2668  {
2669  if (c1 == always || c2 == always)
2670  return always;
2671  if (c1 == nonzero || c2 == nonzero)
2672  return nonzero;
2673  return none;
2674  }
2675 
2676 
2677  // ---------------------- inline and template functions --------------------
2678 
2679  template <int dim, int spacedim>
2680  inline unsigned int
2682  {
2683  return dh.get_fe().dofs_per_cell;
2684  }
2685 
2686 
2687  template <int dim, int spacedim>
2688  inline unsigned int
2690  {
2691  return dh.get_fe().dofs_per_face;
2692  }
2693 
2694 
2695  template <int dim, int spacedim>
2696  inline unsigned int
2698  {
2699  return dh.get_fe().dofs_per_vertex;
2700  }
2701 
2702 
2703  template <int dim, int spacedim>
2704  inline unsigned int
2706  {
2707  return dh.get_fe().n_components();
2708  }
2709 
2710 
2711 
2712  template <int dim, int spacedim>
2713  inline bool
2715  {
2716  return dh.get_fe().is_primitive();
2717  }
2718 
2719 
2720  template <int dim, int spacedim>
2721  inline unsigned int
2723  {
2724  return dh.get_fe_collection().max_dofs_per_cell();
2725  }
2726 
2727 
2728  template <int dim, int spacedim>
2729  inline unsigned int
2731  {
2732  return dh.get_fe_collection().max_dofs_per_face();
2733  }
2734 
2735 
2736  template <int dim, int spacedim>
2737  inline unsigned int
2739  {
2740  return dh.get_fe_collection().max_dofs_per_vertex();
2741  }
2742 
2743 
2744  template <int dim, int spacedim>
2745  inline unsigned int
2747  {
2748  return dh.get_fe(0).n_components();
2749  }
2750 
2751 
2752  template <int dim, int spacedim>
2753  inline bool
2755  {
2756  return dh.get_fe(0).is_primitive();
2757  }
2758 
2759 
2760  template <typename DoFHandlerType, class Comp>
2761  void
2764  & mapping,
2765  const DoFHandlerType &dof_handler,
2768  Comp> & point_to_index_map)
2769  {
2770  // let the checking of arguments be
2771  // done by the function first
2772  // called
2773  std::vector<Point<DoFHandlerType::space_dimension>> support_points(
2774  dof_handler.n_dofs());
2775  map_dofs_to_support_points(mapping, dof_handler, support_points);
2776  // now copy over the results of the
2777  // previous function into the
2778  // output arg
2779  point_to_index_map.clear();
2780  for (types::global_dof_index i = 0; i < dof_handler.n_dofs(); ++i)
2781  point_to_index_map[support_points[i]] = i;
2782  }
2783 } // namespace DoFTools
2784 
2785 #endif
2786 
2787 DEAL_II_NAMESPACE_CLOSE
2788 
2789 #endif
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:2300
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:2495
unsigned int max_dofs_per_vertex(const DoFHandler< dim, spacedim > &dh)
std::vector< IndexSet > locally_owned_dofs_per_subdomain(const DoFHandlerType &dof_handler)
Definition: dof_tools.cc:1331
const types::subdomain_id invalid_subdomain_id
Definition: types.h:258
static ::ExceptionBase & ExcGridsDontMatch()
void make_boundary_sparsity_pattern(const DoFHandlerType &dof, const std::vector< types::global_dof_index > &dof_to_boundary_mapping, SparsityPatternType &sparsity_pattern)
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:2402
void make_periodicity_constraints(const FaceIterator &face_1, const typename identity< FaceIterator >::type &face_2, AffineConstraints< double > &constraints, const ComponentMask &component_mask=ComponentMask(), const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const FullMatrix< double > &matrix=FullMatrix< double >(), const std::vector< unsigned int > &first_vector_components=std::vector< unsigned int >())
void get_subdomain_association(const DoFHandlerType &dof_handler, std::vector< types::subdomain_id > &subdomain)
Definition: dof_tools.cc:1506
std::vector< types::global_dof_index > get_dofs_on_patch(const std::vector< typename DoFHandlerType::active_cell_iterator > &patch)
Definition: dof_tools.cc:2775
const hp::FECollection< dim, spacedim > & get_fe() const
void map_support_points_to_dofs(const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &mapping, const DoFHandlerType &dof_handler, std::map< Point< DoFHandlerType::space_dimension >, types::global_dof_index, Comp > &point_to_index_map)
bool fe_is_primitive(const DoFHandler< dim, spacedim > &dh)
static ::ExceptionBase & ExcFiniteElementsDontMatch()
unsigned int count_dofs_with_subdomain_association(const DoFHandlerType &dof_handler, const types::subdomain_id subdomain)
Definition: dof_tools.cc:1605
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
void map_dof_to_boundary_indices(const DoFHandlerType &dof_handler, std::vector< types::global_dof_index > &mapping)
Definition: dof_tools.cc:2041
Table< 2, Coupling > dof_couplings_from_component_couplings(const FiniteElement< dim, spacedim > &fe, const Table< 2, Coupling > &component_couplings)
void extract_locally_relevant_level_dofs(const DoFHandlerType &dof_handler, const unsigned int level, IndexSet &dof_set)
Definition: dof_tools.cc:1179
void make_hanging_node_constraints(const DoFHandlerType &dof_handler, AffineConstraints< number > &constraints)
LinearAlgebra::distributed::Vector< Number > Vector
IndexSet dof_indices_with_subdomain_association(const DoFHandlerType &dof_handler, const types::subdomain_id subdomain)
Definition: dof_tools.cc:1621
UpdateFlags operator|(const UpdateFlags f1, const UpdateFlags f2)
Definition: tria.h:81
unsigned int max_dofs_per_face(const DoFHandler< dim, spacedim > &dh)
unsigned int subdomain_id
Definition: types.h:43
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:2371
static ::ExceptionBase & ExcGridNotCoarser()
void extract_locally_owned_dofs(const DoFHandlerType &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:1090
Abstract base class for mapping classes.
Definition: dof_tools.h:57
unsigned int max_dofs_per_cell(const DoFHandler< dim, spacedim > &dh)
#define DeclException0(Exception0)
Definition: exceptions.h:473
void map_dofs_to_support_points(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, std::vector< Point< spacedim >> &support_points)
Definition: dof_tools.cc:2222
void make_flux_sparsity_pattern(const DoFHandlerType &dof_handler, SparsityPatternType &sparsity_pattern)
void compute_intergrid_transfer_representation(const DoFHandler< dim, spacedim > &coarse_grid, const unsigned int coarse_component, const DoFHandler< dim, spacedim > &fine_grid, const unsigned int fine_component, const InterGridMap< DoFHandler< dim, spacedim >> &coarse_to_fine_grid_map, std::vector< std::map< types::global_dof_index, float >> &transfer_representation)
void extract_locally_relevant_dofs(const DoFHandlerType &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:1137
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:591
void extract_hanging_node_dofs(const DoFHandler< dim, spacedim > &dof_handler, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:1031
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:730
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:2567
const hp::FECollection< dim, spacedim > & get_fe_collection() const
Definition: hp.h:117
void make_zero_boundary_constraints(const DoFHandlerType< dim, spacedim > &dof, const types::boundary_id boundary_id, AffineConstraints< number > &zero_boundary_constraints, const ComponentMask &component_mask=ComponentMask())
unsigned int global_dof_index
Definition: types.h:89
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
void extract_constant_modes(const DoFHandlerType &dof_handler, const ComponentMask &component_mask, std::vector< std::vector< bool >> &constant_modes)
Definition: dof_tools.cc:1228
std::vector< IndexSet > locally_relevant_dofs_per_subdomain(const DoFHandlerType &dof_handler)
Definition: dof_tools.cc:1431
void extract_locally_active_dofs(const DoFHandlerType &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:1102
void extract_dofs(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:402
void extract_subdomain_dofs(const DoFHandlerType &dof_handler, const types::subdomain_id subdomain_id, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:1057
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:1953
unsigned int count_dofs_on_patch(const std::vector< typename DoFHandlerType::active_cell_iterator > &patch)
Definition: dof_tools.cc:2746
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
Definition: table.h:37
void extract_level_dofs(const unsigned int level, const DoFHandlerType &dof, const ComponentMask &component_mask, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:520
void compute_intergrid_constraints(const DoFHandler< dim, spacedim > &coarse_grid, const unsigned int coarse_component, const DoFHandler< dim, spacedim > &fine_grid, const unsigned int fine_component, const InterGridMap< DoFHandler< dim, spacedim >> &coarse_to_fine_grid_map, AffineConstraints< double > &constraints)
void make_sparsity_pattern(const DoFHandlerType &dof_handler, SparsityPatternType &sparsity_pattern, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
void get_active_fe_indices(const DoFHandlerType &dof_handler, std::vector< unsigned int > &active_fe_indices)
Definition: dof_tools.cc:1316
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:317
unsigned int boundary_id
Definition: types.h:111
static ::ExceptionBase & ExcNoFESelected()
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:817
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:2448
void count_dofs_per_component(const DoFHandlerType &dof_handler, std::vector< types::global_dof_index > &dofs_per_component, const bool vector_valued_once=false, std::vector< unsigned int > target_component={})
Definition: dof_tools.cc:1857