16 #ifndef dealii_dof_tools_h 17 #define dealii_dof_tools_h 20 #include <deal.II/base/config.h> 22 #include <deal.II/base/exceptions.h> 23 #include <deal.II/base/index_set.h> 24 #include <deal.II/base/point.h> 26 #include <deal.II/dofs/dof_handler.h> 28 #include <deal.II/fe/component_mask.h> 30 #include <deal.II/hp/dof_handler.h> 32 #include <deal.II/lac/affine_constraints.h> 40 DEAL_II_NAMESPACE_OPEN
43 template <
int dim,
typename RangeNumberType>
45 template <
int dim,
int spacedim>
49 template <
int dim,
int spacedim>
51 template <
int dim,
int spacedim>
54 template <
class MeshType>
56 template <
int dim,
int spacedim>
59 template <
int dim,
class T>
61 template <
typename Number>
66 template <
typename CellIterator>
67 struct PeriodicFacePair;
237 template <
int dim,
int spacedim>
252 template <
int dim,
int spacedim>
263 template <
int dim,
int spacedim>
276 template <
int dim,
int spacedim>
277 std::vector<Table<2, Coupling>>
299 template <
int dim,
int spacedim>
300 DEAL_II_DEPRECATED
unsigned int 310 template <
int dim,
int spacedim>
311 DEAL_II_DEPRECATED
unsigned int 325 template <
int dim,
int spacedim>
326 DEAL_II_DEPRECATED
unsigned int 339 template <
int dim,
int spacedim>
340 DEAL_II_DEPRECATED
unsigned int 353 template <
int dim,
int spacedim>
354 DEAL_II_DEPRECATED
unsigned int 367 template <
int dim,
int spacedim>
368 DEAL_II_DEPRECATED
unsigned int 382 template <
int dim,
int spacedim>
383 DEAL_II_DEPRECATED
unsigned int 397 template <
int dim,
int spacedim>
398 DEAL_II_DEPRECATED
unsigned int 412 template <
int dim,
int spacedim>
413 DEAL_II_DEPRECATED
bool 427 template <
int dim,
int spacedim>
428 DEAL_II_DEPRECATED
bool 556 template <
typename DoFHandlerType,
557 typename SparsityPatternType,
558 typename number =
double>
561 const DoFHandlerType & dof_handler,
562 SparsityPatternType & sparsity_pattern,
564 const bool keep_constrained_dofs =
true,
632 template <
typename DoFHandlerType,
633 typename SparsityPatternType,
634 typename number =
double>
637 const DoFHandlerType & dof_handler,
639 SparsityPatternType & sparsity_pattern,
641 const bool keep_constrained_dofs =
true,
664 template <
typename DoFHandlerType,
typename SparsityPatternType>
667 const DoFHandlerType &dof_col,
668 SparsityPatternType & sparsity);
715 template <
typename DoFHandlerType,
typename SparsityPatternType>
718 SparsityPatternType & sparsity_pattern);
728 template <
typename DoFHandlerType,
729 typename SparsityPatternType,
733 const DoFHandlerType & dof_handler,
734 SparsityPatternType & sparsity_pattern,
736 const bool keep_constrained_dofs =
true,
759 template <
typename DoFHandlerType,
typename SparsityPatternType>
762 const DoFHandlerType & dof,
763 SparsityPatternType & sparsity,
777 template <
typename DoFHandlerType,
778 typename SparsityPatternType,
782 SparsityPatternType & sparsity,
784 const bool keep_constrained_dofs,
798 template <
typename DoFHandlerType,
typename SparsityPatternType>
801 const DoFHandlerType & dof,
802 const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
803 SparsityPatternType & sparsity_pattern);
822 template <
typename DoFHandlerType,
823 typename SparsityPatternType,
827 const DoFHandlerType &dof,
831 const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
832 SparsityPatternType & sparsity);
884 template <
typename DoFHandlerType,
typename number>
956 template <
int dim,
int spacedim>
960 const unsigned int coarse_component,
962 const unsigned int fine_component,
983 template <
int dim,
int spacedim>
987 const unsigned int coarse_component,
989 const unsigned int fine_component,
991 std::vector<std::map<types::global_dof_index, float>>
992 &transfer_representation);
1156 template <
typename FaceIterator>
1159 const FaceIterator & face_1,
1160 const typename identity<FaceIterator>::type &face_2,
1163 const bool face_orientation =
true,
1164 const bool face_flip =
false,
1165 const bool face_rotation =
false,
1167 const std::vector<unsigned int> &first_vector_components =
1168 std::vector<unsigned int>());
1193 template <
typename DoFHandlerType>
1201 const std::vector<unsigned int> &first_vector_components =
1202 std::vector<unsigned int>());
1236 template <
typename DoFHandlerType>
1239 const DoFHandlerType & dof_handler,
1242 const int direction,
1272 template <
typename DoFHandlerType>
1275 const DoFHandlerType & dof_handler,
1277 const int direction,
1305 template <
int dim,
int spacedim>
1308 std::vector<bool> & selected_dofs);
1314 template <
int dim,
int spacedim>
1351 template <
int dim,
int spacedim>
1355 std::vector<bool> & selected_dofs);
1360 template <
int dim,
int spacedim>
1364 std::vector<bool> & selected_dofs);
1390 template <
int dim,
int spacedim>
1394 std::vector<bool> & selected_dofs);
1399 template <
int dim,
int spacedim>
1403 std::vector<bool> & selected_dofs);
1409 template <
typename DoFHandlerType>
1412 const DoFHandlerType &dof,
1414 std::vector<bool> & selected_dofs);
1420 template <
typename DoFHandlerType>
1423 const DoFHandlerType &dof,
1425 std::vector<bool> & selected_dofs);
1479 template <
typename DoFHandlerType>
1483 std::vector<bool> & selected_dofs,
1484 const std::set<types::boundary_id> &boundary_ids =
1485 std::set<types::boundary_id>());
1517 template <
typename DoFHandlerType>
1522 const std::set<types::boundary_id> &boundary_ids =
1523 std::set<types::boundary_id>());
1541 template <
typename DoFHandlerType>
1544 const DoFHandlerType & dof_handler,
1546 std::vector<bool> & selected_dofs,
1547 const std::set<types::boundary_id> &boundary_ids =
1548 std::set<types::boundary_id>());
1577 template <
typename DoFHandlerType,
typename number =
double>
1580 const DoFHandlerType &dof_handler,
1581 const std::function<
1582 bool(
const typename DoFHandlerType::active_cell_iterator &)> &predicate,
1618 template <
typename DoFHandlerType>
1622 std::vector<std::vector<bool>> &constant_modes);
1638 template <
typename DoFHandlerType>
1642 std::vector<bool> & selected_dofs);
1654 template <
typename DoFHandlerType>
1655 DEAL_II_DEPRECATED
void 1674 template <
typename DoFHandlerType>
1688 template <
typename DoFHandlerType>
1709 template <
typename DoFHandlerType>
1710 std::vector<IndexSet>
1729 template <
typename DoFHandlerType>
1730 std::vector<IndexSet>
1738 template <
typename DoFHandlerType>
1741 const unsigned int level,
1775 template <
typename DoFHandlerType>
1778 std::vector<types::subdomain_id> &subdomain);
1805 template <
typename DoFHandlerType>
1829 template <
typename DoFHandlerType>
1832 const DoFHandlerType & dof_handler,
1834 std::vector<unsigned int> &n_dofs_on_subdomain);
1857 template <
typename DoFHandlerType>
1923 template <
typename DoFHandlerType>
1924 std::vector<types::global_dof_index>
1926 const std::vector<typename DoFHandlerType::active_cell_iterator> &patch);
1948 template <
int dim,
int spacedim>
1952 const unsigned int level,
1953 const std::vector<bool> & selected_dofs = {},
2007 template <
typename DoFHandlerType>
2008 std::vector<unsigned int>
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);
2032 template <
typename DoFHandlerType>
2033 std::vector<unsigned int>
2035 const DoFHandlerType &dof_handler,
2036 const unsigned int level,
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);
2080 template <
typename DoFHandlerType>
2083 const DoFHandlerType &dof_handler,
2084 const unsigned int level,
2085 const bool interior_dofs_only,
2086 const bool boundary_dofs =
false);
2107 template <
typename DoFHandlerType>
2110 const DoFHandlerType &dof_handler,
2111 const unsigned int level,
2112 const bool interior_dofs_only =
false);
2156 template <
typename DoFHandlerType>
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 = {});
2180 template <
typename DoFHandlerType>
2183 std::vector<types::global_dof_index> &dofs_per_block,
2184 const std::vector<unsigned int> & target_block =
2185 std::vector<unsigned int>());
2197 template <
typename DoFHandlerType>
2200 std::vector<unsigned int> &active_fe_indices);
2236 template <
typename DoFHandlerType>
2239 const std::vector<typename DoFHandlerType::active_cell_iterator> &patch);
2263 template <
typename DoFHandlerType>
2266 std::vector<types::global_dof_index> &mapping);
2278 template <
typename DoFHandlerType>
2281 const std::set<types::boundary_id> &boundary_ids,
2282 std::vector<types::global_dof_index> &mapping);
2301 template <
int dim,
int spacedim>
2311 template <
int dim,
int spacedim>
2314 const ::hp::MappingCollection<dim, spacedim> &mapping,
2345 template <
int dim,
int spacedim>
2355 template <
int dim,
int spacedim>
2358 const ::hp::MappingCollection<dim, spacedim> &mapping,
2383 template <
typename DoFHandlerType,
class Comp>
2388 const DoFHandlerType &dof_handler,
2391 Comp> & point_to_index_map);
2429 template <
typename DoFHandlerType,
typename Number>
2432 const Vector<Number> &cell_data,
2434 const unsigned int component = 0);
2513 template <
int spacedim>
2562 template <
int,
int>
class DoFHandlerType,
2566 const DoFHandlerType<dim, spacedim> &dof,
2583 template <
int,
int>
class DoFHandlerType,
2587 const DoFHandlerType<dim, spacedim> &dof,
2679 template <
int dim,
int spacedim>
2683 return dh.
get_fe().dofs_per_cell;
2687 template <
int dim,
int spacedim>
2691 return dh.
get_fe().dofs_per_face;
2695 template <
int dim,
int spacedim>
2699 return dh.
get_fe().dofs_per_vertex;
2703 template <
int dim,
int spacedim>
2707 return dh.
get_fe().n_components();
2712 template <
int dim,
int spacedim>
2716 return dh.
get_fe().is_primitive();
2720 template <
int dim,
int spacedim>
2728 template <
int dim,
int spacedim>
2736 template <
int dim,
int spacedim>
2744 template <
int dim,
int spacedim>
2748 return dh.
get_fe(0).n_components();
2752 template <
int dim,
int spacedim>
2756 return dh.
get_fe(0).is_primitive();
2760 template <
typename DoFHandlerType,
class Comp>
2765 const DoFHandlerType &dof_handler,
2768 Comp> & point_to_index_map)
2773 std::vector<Point<DoFHandlerType::space_dimension>> support_points(
2774 dof_handler.n_dofs());
2779 point_to_index_map.clear();
2781 point_to_index_map[support_points[i]] = i;
2787 DEAL_II_NAMESPACE_CLOSE
const types::subdomain_id invalid_subdomain_id
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)
const hp::FECollection< dim, spacedim > & get_fe() const
static ::ExceptionBase & ExcFiniteElementsDontMatch()
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
void make_hanging_node_constraints(const DoFHandlerType &dof_handler, AffineConstraints< number > &constraints)
LinearAlgebra::distributed::Vector< Number > Vector
UpdateFlags operator|(const UpdateFlags f1, const UpdateFlags f2)
unsigned int subdomain_id
static ::ExceptionBase & ExcGridNotCoarser()
Abstract base class for mapping classes.
#define DeclException0(Exception0)
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)
const hp::FECollection< dim, spacedim > & get_fe_collection() const
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
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
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)
static ::ExceptionBase & ExcNoFESelected()