16 #ifndef dealii_dof_tools_h 17 #define dealii_dof_tools_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/exceptions.h> 22 #include <deal.II/base/index_set.h> 23 #include <deal.II/base/point.h> 24 #include <deal.II/lac/constraint_matrix.h> 25 #include <deal.II/dofs/dof_handler.h> 26 #include <deal.II/fe/component_mask.h> 27 #include <deal.II/hp/dof_handler.h> 35 DEAL_II_NAMESPACE_OPEN
38 template <
int dim,
typename RangeNumberType>
class Function;
46 template <
int dim,
int spacedim>
class Mapping;
48 template <
int dim,
class T>
class Table;
49 template <
typename Number>
class Vector;
53 template <
typename CellIterator>
struct PeriodicFacePair;
222 template <
int dim,
int spacedim>
237 template <
int dim,
int spacedim>
248 template <
int dim,
int spacedim>
260 template <
int dim,
int spacedim>
261 std::vector<Table<2,Coupling> >
279 template <
int dim,
int spacedim>
291 template <
int dim,
int spacedim>
307 template <
int dim,
int spacedim>
322 template <
int dim,
int spacedim>
337 template <
int dim,
int spacedim>
352 template <
int dim,
int spacedim>
368 template <
int dim,
int spacedim>
384 template <
int dim,
int spacedim>
400 template <
int dim,
int spacedim>
416 template <
int dim,
int spacedim>
545 template <
typename DoFHandlerType,
typename SparsityPatternType>
548 SparsityPatternType &sparsity_pattern,
550 const bool keep_constrained_dofs =
true,
618 template <
typename DoFHandlerType,
typename SparsityPatternType>
622 SparsityPatternType &sparsity_pattern,
624 const bool keep_constrained_dofs =
true,
647 template <
typename DoFHandlerType,
typename SparsityPatternType>
650 const DoFHandlerType &dof_col,
651 SparsityPatternType &sparsity);
698 template <
typename DoFHandlerType,
typename SparsityPatternType>
701 SparsityPatternType &sparsity_pattern);
711 template <
typename DoFHandlerType,
typename SparsityPatternType>
714 SparsityPatternType &sparsity_pattern,
716 const bool keep_constrained_dofs =
true,
738 template <
typename DoFHandlerType,
typename SparsityPatternType>
741 SparsityPatternType &sparsity,
753 template <
typename DoFHandlerType,
typename SparsityPatternType>
756 SparsityPatternType &sparsity,
758 const bool keep_constrained_dofs,
772 template <
typename DoFHandlerType,
typename SparsityPatternType>
775 const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
776 SparsityPatternType &sparsity_pattern);
794 template <
typename DoFHandlerType,
typename SparsityPatternType,
typename number>
797 (
const DoFHandlerType &dof,
799 const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
800 SparsityPatternType &sparsity);
852 template <
typename DoFHandlerType>
924 template <
int dim,
int spacedim>
927 const unsigned int coarse_component,
929 const unsigned int fine_component,
950 template <
int dim,
int spacedim>
953 const unsigned int coarse_component,
955 const unsigned int fine_component,
957 std::vector<std::map<types::global_dof_index, float> > &transfer_representation);
1121 template <
typename FaceIterator>
1124 (
const FaceIterator &face_1,
1125 const typename identity<FaceIterator>::type &face_2,
1128 const bool face_orientation =
true,
1129 const bool face_flip =
false,
1130 const bool face_rotation =
false,
1132 const std::vector<unsigned int> &first_vector_components = std::vector<unsigned int>());
1157 template <
typename DoFHandlerType>
1164 const std::vector<unsigned int> &first_vector_components = std::vector<unsigned int>());
1198 template <
typename DoFHandlerType>
1201 (
const DoFHandlerType &dof_handler,
1204 const int direction,
1234 template <
typename DoFHandlerType>
1237 (
const DoFHandlerType &dof_handler,
1239 const int direction,
1267 template <
int dim,
int spacedim>
1270 std::vector<bool> &selected_dofs);
1276 template <
int dim,
int spacedim>
1313 template <
int dim,
int spacedim>
1317 std::vector<bool> &selected_dofs);
1322 template <
int dim,
int spacedim>
1326 std::vector<bool> &selected_dofs);
1352 template <
int dim,
int spacedim>
1356 std::vector<bool> &selected_dofs);
1361 template <
int dim,
int spacedim>
1365 std::vector<bool> &selected_dofs);
1371 template <
typename DoFHandlerType>
1374 const DoFHandlerType &dof,
1376 std::vector<bool> &selected_dofs);
1382 template <
typename DoFHandlerType>
1385 const DoFHandlerType &dof,
1387 std::vector<bool> &selected_dofs);
1440 template <
typename DoFHandlerType>
1444 std::vector<bool> &selected_dofs,
1445 const std::set<types::boundary_id> &boundary_ids = std::set<types::boundary_id>());
1477 template <
typename DoFHandlerType>
1482 const std::set<types::boundary_id> &boundary_ids = std::set<types::boundary_id>());
1500 template <
typename DoFHandlerType>
1504 std::vector<bool> &selected_dofs,
1505 const std::set<types::boundary_id> &boundary_ids = std::set<types::boundary_id>());
1534 template <
typename DoFHandlerType>
1537 const std::function<
bool(
const typename DoFHandlerType::active_cell_iterator &)> &predicate,
1573 template <
typename DoFHandlerType>
1577 std::vector<std::vector<bool> > &constant_modes);
1593 template <
typename DoFHandlerType>
1597 std::vector<bool> &selected_dofs);
1609 template <
typename DoFHandlerType>
1630 template <
typename DoFHandlerType>
1644 template <
typename DoFHandlerType>
1665 template <
typename DoFHandlerType>
1666 std::vector<IndexSet>
1685 template <
typename DoFHandlerType>
1686 std::vector<IndexSet>
1694 template <
typename DoFHandlerType>
1697 const unsigned int level,
1731 template <
typename DoFHandlerType>
1734 std::vector<types::subdomain_id> &subdomain);
1761 template <
typename DoFHandlerType>
1785 template <
typename DoFHandlerType>
1789 std::vector<unsigned int> &n_dofs_on_subdomain);
1812 template <
typename DoFHandlerType>
1878 template <
typename DoFHandlerType>
1879 std::vector<types::global_dof_index>
1880 get_dofs_on_patch (
const std::vector<typename DoFHandlerType::active_cell_iterator> &patch);
1902 template <
int dim,
int spacedim>
1905 const unsigned int level,
1906 const std::vector<bool> &selected_dofs = std::vector<bool>(),
1959 template <
typename DoFHandlerType>
1960 std::vector<unsigned int>
1962 const DoFHandlerType &dof_handler,
1963 const unsigned int level,
1964 const bool interior_dofs_only,
1965 const bool boundary_patches =
false,
1966 const bool level_boundary_patches =
false,
1967 const bool single_cell_patches =
false,
1968 const bool invert_vertex_mapping =
false);
1983 template <
typename DoFHandlerType>
1984 std::vector<unsigned int>
1986 const DoFHandlerType &dof_handler,
1987 const unsigned int level,
1989 const bool boundary_patches =
false,
1990 const bool level_boundary_patches =
false,
1991 const bool single_cell_patches =
false,
1992 const bool invert_vertex_mapping =
false);
2031 template <
typename DoFHandlerType>
2033 const DoFHandlerType &dof_handler,
2034 const unsigned int level,
2035 const bool interior_dofs_only,
2036 const bool boundary_dofs =
false);
2057 template <
typename DoFHandlerType>
2059 const DoFHandlerType &dof_handler,
2060 const unsigned int level,
2061 const bool interior_dofs_only =
false);
2105 template <
typename DoFHandlerType>
2108 std::vector<types::global_dof_index> &dofs_per_component,
2109 const bool vector_valued_once =
false,
2110 std::vector<unsigned int> target_component
2111 = std::vector<unsigned int>());
2129 template <
typename DoFHandlerType>
2132 std::vector<types::global_dof_index> &dofs_per_block,
2133 const std::vector<unsigned int> &target_block
2134 = std::vector<unsigned int>());
2146 template <
typename DoFHandlerType>
2149 std::vector<unsigned int> &active_fe_indices);
2185 template <
typename DoFHandlerType>
2187 count_dofs_on_patch (
const std::vector<typename DoFHandlerType::active_cell_iterator> &patch);
2211 template <
typename DoFHandlerType>
2214 std::vector<types::global_dof_index> &mapping);
2226 template <
typename DoFHandlerType>
2229 const std::set<types::boundary_id> &boundary_ids,
2230 std::vector<types::global_dof_index> &mapping);
2249 template <
int dim,
int spacedim>
2259 template <
int dim,
int spacedim>
2292 template <
int dim,
int spacedim>
2301 template <
int dim,
int spacedim>
2328 template <
typename DoFHandlerType,
class Comp>
2332 const DoFHandlerType &dof_handler,
2371 template <
typename DoFHandlerType,
typename Number>
2374 const Vector<Number> &cell_data,
2376 const unsigned int component = 0);
2451 template <
int spacedim>
2497 template <
int dim,
int spacedim,
template <
int,
int>
class DoFHandlerType>
2514 template <
int dim,
int spacedim,
template <
int,
int>
class DoFHandlerType>
2610 template <
int dim,
int spacedim>
2615 return dh.
get_fe().dofs_per_cell;
2619 template <
int dim,
int spacedim>
2624 return dh.
get_fe().dofs_per_face;
2628 template <
int dim,
int spacedim>
2633 return dh.
get_fe().dofs_per_vertex;
2637 template <
int dim,
int spacedim>
2642 return dh.
get_fe().n_components();
2647 template <
int dim,
int spacedim>
2652 return dh.
get_fe().is_primitive();
2656 template <
int dim,
int spacedim>
2665 template <
int dim,
int spacedim>
2674 template <
int dim,
int spacedim>
2683 template <
int dim,
int spacedim>
2688 return dh.
get_fe(0).n_components();
2692 template <
int dim,
int spacedim>
2697 return dh.
get_fe(0).is_primitive();
2701 template <
typename DoFHandlerType,
class Comp>
2706 const DoFHandlerType &dof_handler,
2712 std::vector<Point<DoFHandlerType::space_dimension> > support_points (dof_handler.n_dofs());
2717 point_to_index_map.clear ();
2719 point_to_index_map[support_points[i]] = i;
2725 DEAL_II_NAMESPACE_CLOSE
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 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)
static ::ExceptionBase & ExcFiniteElementsDontMatch()
UpdateFlags operator|(const UpdateFlags f1, const UpdateFlags f2)
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, ConstraintMatrix &constraints)
const hp::FECollection< dim, spacedim > & get_fe() const
void make_hanging_node_constraints(const DoFHandlerType &dof_handler, ConstraintMatrix &constraints)
static ::ExceptionBase & ExcGridNotCoarser()
unsigned int global_dof_index
Abstract base class for mapping classes.
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
#define DeclException0(Exception0)
void make_flux_sparsity_pattern(const DoFHandlerType &dof_handler, SparsityPatternType &sparsity_pattern)
unsigned int subdomain_id
const hp::FECollection< dim, spacedim > & get_fe_collection() const
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
void make_sparsity_pattern(const DoFHandlerType &dof_handler, SparsityPatternType &sparsity_pattern, const ConstraintMatrix &constraints=ConstraintMatrix(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
void make_zero_boundary_constraints(const DoFHandlerType< dim, spacedim > &dof, const types::boundary_id boundary_id, ConstraintMatrix &zero_boundary_constraints, const ComponentMask &component_mask=ComponentMask())
static ::ExceptionBase & ExcNoFESelected()