Reference documentation for deal.II version 9.0.0
|
Enumerations | |
enum | Coupling { none, always, nonzero } |
Functions | |
DoF couplings | |
template<int dim, int spacedim> | |
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) |
template<int dim, int spacedim> | |
void | convert_couplings_to_blocks (const DoFHandler< dim, spacedim > &dof_handler, const Table< 2, Coupling > &table_by_component, std::vector< Table< 2, Coupling > > &tables_by_block) |
template<int dim, int spacedim> | |
Table< 2, Coupling > | dof_couplings_from_component_couplings (const FiniteElement< dim, spacedim > &fe, const Table< 2, Coupling > &component_couplings) |
template<int dim, int spacedim> | |
std::vector< Table< 2, Coupling > > | dof_couplings_from_component_couplings (const hp::FECollection< dim, spacedim > &fe, const Table< 2, Coupling > &component_couplings) |
Functions to support code that generically uses both DoFHandler and hp::DoFHandler | |
template<int dim, int spacedim> | |
unsigned int | max_dofs_per_cell (const DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
unsigned int | max_dofs_per_cell (const hp::DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
unsigned int | max_dofs_per_face (const DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
unsigned int | max_dofs_per_face (const hp::DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
unsigned int | max_dofs_per_vertex (const DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
unsigned int | max_dofs_per_vertex (const hp::DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
unsigned int | n_components (const DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
unsigned int | n_components (const hp::DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
bool | fe_is_primitive (const DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
bool | fe_is_primitive (const hp::DoFHandler< dim, spacedim > &dh) |
Sparsity pattern generation | |
template<typename DoFHandlerType , typename SparsityPatternType > | |
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) |
template<typename DoFHandlerType , typename SparsityPatternType > | |
void | make_sparsity_pattern (const DoFHandlerType &dof_handler, const Table< 2, Coupling > &coupling, SparsityPatternType &sparsity_pattern, const ConstraintMatrix &constraints=ConstraintMatrix(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id) |
template<typename DoFHandlerType , typename SparsityPatternType > | |
void | make_sparsity_pattern (const DoFHandlerType &dof_row, const DoFHandlerType &dof_col, SparsityPatternType &sparsity) |
template<typename DoFHandlerType , typename SparsityPatternType > | |
void | make_flux_sparsity_pattern (const DoFHandlerType &dof_handler, SparsityPatternType &sparsity_pattern) |
template<typename DoFHandlerType , typename SparsityPatternType > | |
void | make_flux_sparsity_pattern (const DoFHandlerType &dof_handler, SparsityPatternType &sparsity_pattern, const ConstraintMatrix &constraints, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id) |
template<typename DoFHandlerType , typename SparsityPatternType > | |
void | make_flux_sparsity_pattern (const DoFHandlerType &dof, SparsityPatternType &sparsity, const Table< 2, Coupling > &cell_integrals_mask, const Table< 2, Coupling > &face_integrals_mask, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id) |
template<typename DoFHandlerType , typename SparsityPatternType > | |
void | make_flux_sparsity_pattern (const DoFHandlerType &dof, SparsityPatternType &sparsity, const ConstraintMatrix &constraints, const bool keep_constrained_dofs, const Table< 2, Coupling > &couplings, const Table< 2, Coupling > &face_couplings, const types::subdomain_id subdomain_id) |
template<typename DoFHandlerType , typename SparsityPatternType > | |
void | make_boundary_sparsity_pattern (const DoFHandlerType &dof, const std::vector< types::global_dof_index > &dof_to_boundary_mapping, SparsityPatternType &sparsity_pattern) |
template<typename DoFHandlerType , typename SparsityPatternType , typename number > | |
void | make_boundary_sparsity_pattern (const DoFHandlerType &dof, const std::map< types::boundary_id, const Function< DoFHandlerType::space_dimension, number > *> &boundary_ids, const std::vector< types::global_dof_index > &dof_to_boundary_mapping, SparsityPatternType &sparsity) |
Hanging nodes and other constraints | |
template<typename DoFHandlerType > | |
void | make_hanging_node_constraints (const DoFHandlerType &dof_handler, ConstraintMatrix &constraints) |
template<int dim, int spacedim> | |
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) |
template<int dim, int spacedim> | |
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) |
Periodic boundary conditions | |
template<typename FaceIterator > | |
void | make_periodicity_constraints (const FaceIterator &face_1, const typename identity< FaceIterator >::type &face_2, ::ConstraintMatrix &constraint_matrix, 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 >()) |
template<typename DoFHandlerType > | |
void | make_periodicity_constraints (const std::vector< GridTools::PeriodicFacePair< typename DoFHandlerType::cell_iterator > > &periodic_faces, ::ConstraintMatrix &constraint_matrix, const ComponentMask &component_mask=ComponentMask(), const std::vector< unsigned int > &first_vector_components=std::vector< unsigned int >()) |
template<typename DoFHandlerType > | |
void | make_periodicity_constraints (const DoFHandlerType &dof_handler, const types::boundary_id b_id1, const types::boundary_id b_id2, const int direction, ::ConstraintMatrix &constraint_matrix, const ComponentMask &component_mask=ComponentMask()) |
template<typename DoFHandlerType > | |
void | make_periodicity_constraints (const DoFHandlerType &dof_handler, const types::boundary_id b_id, const int direction, ::ConstraintMatrix &constraint_matrix, const ComponentMask &component_mask=ComponentMask()) |
Identifying subsets of degrees of freedom with particular properties | |
template<int dim, int spacedim> | |
void | extract_hanging_node_dofs (const DoFHandler< dim, spacedim > &dof_handler, std::vector< bool > &selected_dofs) |
template<int dim, int spacedim> | |
IndexSet | extract_hanging_node_dofs (const DoFHandler< dim, spacedim > &dof_handler) |
template<int dim, int spacedim> | |
void | extract_dofs (const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< bool > &selected_dofs) |
template<int dim, int spacedim> | |
void | extract_dofs (const hp::DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< bool > &selected_dofs) |
template<int dim, int spacedim> | |
void | extract_dofs (const DoFHandler< dim, spacedim > &dof_handler, const BlockMask &block_mask, std::vector< bool > &selected_dofs) |
template<int dim, int spacedim> | |
void | extract_dofs (const hp::DoFHandler< dim, spacedim > &dof_handler, const BlockMask &block_mask, std::vector< bool > &selected_dofs) |
template<typename DoFHandlerType > | |
void | extract_level_dofs (const unsigned int level, const DoFHandlerType &dof, const ComponentMask &component_mask, std::vector< bool > &selected_dofs) |
template<typename DoFHandlerType > | |
void | extract_level_dofs (const unsigned int level, const DoFHandlerType &dof, const BlockMask &component_mask, std::vector< bool > &selected_dofs) |
template<typename DoFHandlerType > | |
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 >()) |
template<typename DoFHandlerType > | |
void | extract_boundary_dofs (const DoFHandlerType &dof_handler, const ComponentMask &component_mask, IndexSet &selected_dofs, const std::set< types::boundary_id > &boundary_ids=std::set< types::boundary_id >()) |
template<typename DoFHandlerType > | |
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 >()) |
template<typename DoFHandlerType > | |
IndexSet | extract_dofs_with_support_contained_within (const DoFHandlerType &dof_handler, const std::function< bool(const typename DoFHandlerType::active_cell_iterator &)> &predicate, const ConstraintMatrix &cm=ConstraintMatrix()) |
template<typename DoFHandlerType > | |
void | extract_constant_modes (const DoFHandlerType &dof_handler, const ComponentMask &component_mask, std::vector< std::vector< bool > > &constant_modes) |
Parallelization and domain decomposition | |
template<typename DoFHandlerType > | |
void | extract_subdomain_dofs (const DoFHandlerType &dof_handler, const types::subdomain_id subdomain_id, std::vector< bool > &selected_dofs) |
template<typename DoFHandlerType > | |
void | extract_locally_owned_dofs (const DoFHandlerType &dof_handler, IndexSet &dof_set) |
template<typename DoFHandlerType > | |
void | extract_locally_active_dofs (const DoFHandlerType &dof_handler, IndexSet &dof_set) |
template<typename DoFHandlerType > | |
void | extract_locally_relevant_dofs (const DoFHandlerType &dof_handler, IndexSet &dof_set) |
template<typename DoFHandlerType > | |
std::vector< IndexSet > | locally_owned_dofs_per_subdomain (const DoFHandlerType &dof_handler) |
template<typename DoFHandlerType > | |
std::vector< IndexSet > | locally_relevant_dofs_per_subdomain (const DoFHandlerType &dof_handler) |
template<typename DoFHandlerType > | |
void | extract_locally_relevant_level_dofs (const DoFHandlerType &dof_handler, const unsigned int level, IndexSet &dof_set) |
template<typename DoFHandlerType > | |
void | get_subdomain_association (const DoFHandlerType &dof_handler, std::vector< types::subdomain_id > &subdomain) |
template<typename DoFHandlerType > | |
unsigned int | count_dofs_with_subdomain_association (const DoFHandlerType &dof_handler, const types::subdomain_id subdomain) |
template<typename DoFHandlerType > | |
void | count_dofs_with_subdomain_association (const DoFHandlerType &dof_handler, const types::subdomain_id subdomain, std::vector< unsigned int > &n_dofs_on_subdomain) |
template<typename DoFHandlerType > | |
IndexSet | dof_indices_with_subdomain_association (const DoFHandlerType &dof_handler, const types::subdomain_id subdomain) |
DoF indices on patches of cells | |
Create structures containing a large set of degrees of freedom for small patches of cells. The resulting objects can be used in RelaxationBlockSOR and related classes to implement Schwarz preconditioners and smoothers, where the subdomains consist of small numbers of cells only. | |
template<typename DoFHandlerType > | |
std::vector< types::global_dof_index > | get_dofs_on_patch (const std::vector< typename DoFHandlerType::active_cell_iterator > &patch) |
template<int dim, int spacedim> | |
void | make_cell_patches (SparsityPattern &block_list, const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, const std::vector< bool > &selected_dofs=std::vector< bool >(), const types::global_dof_index offset=0) |
template<typename DoFHandlerType > | |
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) |
template<typename DoFHandlerType > | |
std::vector< unsigned int > | make_vertex_patches (SparsityPattern &block_list, const DoFHandlerType &dof_handler, const unsigned int level, const BlockMask &exclude_boundary_dofs=BlockMask(), const bool boundary_patches=false, const bool level_boundary_patches=false, const bool single_cell_patches=false, const bool invert_vertex_mapping=false) |
template<typename DoFHandlerType > | |
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) |
template<typename DoFHandlerType > | |
void | make_single_patch (SparsityPattern &block_list, const DoFHandlerType &dof_handler, const unsigned int level, const bool interior_dofs_only=false) |
Counting degrees of freedom and related functions | |
template<typename DoFHandlerType > | |
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=std::vector< unsigned int >()) |
template<typename DoFHandlerType > | |
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 >()) |
template<typename DoFHandlerType > | |
void | get_active_fe_indices (const DoFHandlerType &dof_handler, std::vector< unsigned int > &active_fe_indices) |
template<typename DoFHandlerType > | |
unsigned int | count_dofs_on_patch (const std::vector< typename DoFHandlerType::active_cell_iterator > &patch) |
Functions that return different DoF mappings | |
template<typename DoFHandlerType > | |
void | map_dof_to_boundary_indices (const DoFHandlerType &dof_handler, std::vector< types::global_dof_index > &mapping) |
template<typename DoFHandlerType > | |
void | map_dof_to_boundary_indices (const DoFHandlerType &dof_handler, const std::set< types::boundary_id > &boundary_ids, std::vector< types::global_dof_index > &mapping) |
template<int dim, int spacedim> | |
void | map_dofs_to_support_points (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, std::vector< Point< spacedim > > &support_points) |
template<int dim, int spacedim> | |
void | map_dofs_to_support_points (const ::hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof_handler, std::vector< Point< spacedim > > &support_points) |
template<int dim, int spacedim> | |
void | map_dofs_to_support_points (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, std::map< types::global_dof_index, Point< spacedim > > &support_points) |
template<int dim, int spacedim> | |
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) |
template<typename DoFHandlerType , class Comp > | |
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) |
Miscellaneous | |
template<typename DoFHandlerType , typename Number > | |
void | distribute_cell_to_dof_vector (const DoFHandlerType &dof_handler, const Vector< Number > &cell_data, Vector< double > &dof_data, const unsigned int component=0) |
template<int spacedim> | |
void | write_gnuplot_dof_support_point_info (std::ostream &out, const std::map< types::global_dof_index, Point< spacedim > > &support_points) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType> | |
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()) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType> | |
void | make_zero_boundary_constraints (const DoFHandlerType< dim, spacedim > &dof, ConstraintMatrix &zero_boundary_constraints, const ComponentMask &component_mask=ComponentMask()) |
Exceptions | |
static ::ExceptionBase & | ExcFiniteElementsDontMatch () |
static ::ExceptionBase & | ExcGridNotCoarser () |
static ::ExceptionBase & | ExcGridsDontMatch () |
static ::ExceptionBase & | ExcNoFESelected () |
static ::ExceptionBase & | ExcInvalidBoundaryIndicator () |
This is a collection of functions operating on, and manipulating the numbers of degrees of freedom. The documentation of the member functions will provide more information, but for functions that exist in multiple versions, there are sections in this global documentation stating some commonalities.
When assembling system matrices, the entries are usually of the form \(a_{ij} = a(\phi_i, \phi_j)\), where \(a\) is a bilinear functional, often an integral. When using sparse matrices, we therefore only need to reserve space for those \(a_{ij}\) only, which are nonzero, which is the same as to say that the basis functions \(\phi_i\) and \(\phi_j\) have a nonempty intersection of their support. Since the support of basis functions is bound only on cells on which they are located or to which they are adjacent, to determine the sparsity pattern it is sufficient to loop over all cells and connect all basis functions on each cell with all other basis functions on that cell. There may be finite elements for which not all basis functions on a cell connect with each other, but no use of this case is made since no examples where this occurs are known to the author.
When projecting the traces of functions to the boundary or parts thereof, one needs to build matrices and vectors that act only on those degrees of freedom that are located on the boundary, rather than on all degrees of freedom. One could do that by simply building matrices in which the entries for all interior DoFs are zero, but such matrices are always very rank deficient and not very practical to work with.
What is needed instead in this case is a numbering of the boundary degrees of freedom, i.e. we should enumerate all the degrees of freedom that are sitting on the boundary, and exclude all other (interior) degrees of freedom. The map_dof_to_boundary_indices() function does exactly this: it provides a vector with as many entries as there are degrees of freedom on the whole domain, with each entry being the number in the numbering of the boundary or numbers::invalid_dof_index if the dof is not on the boundary.
With this vector, one can get, for any given degree of freedom, a unique number among those DoFs that sit on the boundary; or, if your DoF was interior to the domain, the result would be numbers::invalid_dof_index. We need this mapping, for example, to build the mass matrix on the boundary (for this, see make_boundary_sparsity_pattern() function, the corresponding section below, as well as the MatrixCreator namespace documentation).
Actually, there are two map_dof_to_boundary_indices() functions, one producing a numbering for all boundary degrees of freedom and one producing a numbering for only parts of the boundary, namely those parts for which the boundary indicator is listed in a set of indicators given to the function. The latter case is needed if, for example, we would only want to project the boundary values for the Dirichlet part of the boundary. You then give the function a list of boundary indicators referring to Dirichlet parts on which the projection is to be performed. The parts of the boundary on which you want to project need not be contiguous; however, it is not guaranteed that the indices of each of the boundary parts are continuous, i.e. the indices of degrees of freedom on different parts may be intermixed.
Degrees of freedom on the boundary but not on one of the specified boundary parts are given the index numbers::invalid_dof_index, as if they were in the interior. If no boundary indicator was given or if no face of a cell has a boundary indicator contained in the given list, the vector of new indices consists solely of numbers::invalid_dof_index.
(As a side note, for corner cases: The question what a degree of freedom on the boundary is, is not so easy. It should really be a degree of freedom of which the respective basis function has nonzero values on the boundary. At least for Lagrange elements this definition is equal to the statement that the off-point, or what deal.II calls support_point, of the shape function, i.e. the point where the function assumes its nominal value (for Lagrange elements this is the point where it has the function value 1), is located on the boundary. We do not check this directly, the criterion is rather defined through the information the finite element class gives: the FiniteElement class defines the numbers of basis functions per vertex, per line, and so on and the basis functions are numbered after this information; a basis function is to be considered to be on the face of a cell (and thus on the boundary if the cell is at the boundary) according to it belonging to a vertex, line, etc but not to the interior of the cell. The finite element uses the same cell-wise numbering so that we can say that if a degree of freedom was numbered as one of the dofs on lines, we assume that it is located on the line. Where the off-point actually is, is a secret of the finite element (well, you can ask it, but we don't do it here) and not relevant in this context.)
In some cases, one wants to only work with DoFs that sit on the boundary. One application is, for example, if rather than interpolating non- homogenous boundary values, one would like to project them. For this, we need two things: a way to identify nodes that are located on (parts of) the boundary, and a way to build matrices out of only degrees of freedom that are on the boundary (i.e. much smaller matrices, in which we do not even build the large zero block that stems from the fact that most degrees of freedom have no support on the boundary of the domain). The first of these tasks is done by the map_dof_to_boundary_indices() function (described above).
The second part requires us first to build a sparsity pattern for the couplings between boundary nodes, and then to actually build the components of this matrix. While actually computing the entries of these small boundary matrices is discussed in the MatrixCreator namespace, the creation of the sparsity pattern is done by the create_boundary_sparsity_pattern() function. For its work, it needs to have a numbering of all those degrees of freedom that are on those parts of the boundary that we are interested in. You can get this from the map_dof_to_boundary_indices() function. It then builds the sparsity pattern corresponding to integrals like \(\int_\Gamma \varphi_{b2d(i)} \varphi_{b2d(j)} dx\), where \(i\) and \(j\) are indices into the matrix, and \(b2d(i)\) is the global DoF number of a degree of freedom sitting on a boundary (i.e., \(b2d\) is the inverse of the mapping returned by map_dof_to_boundary_indices() function).
enum DoFTools::Coupling |
The flags used in tables by certain make_*_pattern
functions to describe whether two components of the solution couple in the bilinear forms corresponding to cell or face terms. An example of using these flags is shown in the introduction of step-46.
In the descriptions of the individual elements below, remember that these flags are used as elements of tables of size FiniteElement::n_components times FiniteElement::n_components where each element indicates whether two components do or do not couple.
Enumerator | |
---|---|
none | Two components do not couple. |
always | Two components do couple. |
nonzero | Two components couple only if their shape functions are both nonzero on a given face. This flag is only used when computing integrals over faces of cells, e.g., in DoFTools::make_flux_sparsity_pattern(). |
Definition at line 190 of file dof_tools.h.
void DoFTools::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 | ||
) |
Map a coupling table from the user friendly organization by components to the organization by blocks. Specializations of this function for DoFHandler and hp::DoFHandler are required due to the different results of their finite element access.
The return vector will be initialized to the correct length inside this function.
Definition at line 2227 of file dof_tools.cc.
void DoFTools::convert_couplings_to_blocks | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
const Table< 2, Coupling > & | table_by_component, | ||
std::vector< Table< 2, Coupling > > & | tables_by_block | ||
) |
Map a coupling table from the user friendly organization by components to the organization by blocks. Specializations of this function for DoFHandler and hp::DoFHandler are required due to the different results of their finite element access.
The return vector will be initialized to the correct length inside this function.
Definition at line 2202 of file dof_tools.cc.
Table< 2, Coupling > DoFTools::dof_couplings_from_component_couplings | ( | const FiniteElement< dim, spacedim > & | fe, |
const Table< 2, Coupling > & | component_couplings | ||
) |
Given a finite element and a table how the vector components of it couple with each other, compute and return a table that describes how the individual shape functions couple with each other.
Definition at line 624 of file dof_tools_sparsity.cc.
std::vector< Table< 2, Coupling > > DoFTools::dof_couplings_from_component_couplings | ( | const hp::FECollection< dim, spacedim > & | fe, |
const Table< 2, Coupling > & | component_couplings | ||
) |
Same function as above for a collection of finite elements, returning a collection of tables.
The function currently treats DoFTools::Couplings::nonzero the same as DoFTools::Couplings::always .
Definition at line 669 of file dof_tools_sparsity.cc.
unsigned int DoFTools::max_dofs_per_cell | ( | const DoFHandler< dim, spacedim > & | dh | ) |
Maximal number of degrees of freedom on a cell.
dh.get_fe_collection().max_dofs_per_cell()
.unsigned int DoFTools::max_dofs_per_cell | ( | const hp::DoFHandler< dim, spacedim > & | dh | ) |
Maximal number of degrees of freedom on a cell.
dh.get_fe_collection().max_dofs_per_cell()
.unsigned int DoFTools::max_dofs_per_face | ( | const DoFHandler< dim, spacedim > & | dh | ) |
Maximal number of degrees of freedom on a face.
This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.
dh.get_fe_collection().max_dofs_per_face()
.unsigned int DoFTools::max_dofs_per_face | ( | const hp::DoFHandler< dim, spacedim > & | dh | ) |
Maximal number of degrees of freedom on a face.
This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.
dh.get_fe_collection().max_dofs_per_face()
.unsigned int DoFTools::max_dofs_per_vertex | ( | const DoFHandler< dim, spacedim > & | dh | ) |
Maximal number of degrees of freedom on a vertex.
This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.
dh.get_fe_collection().max_dofs_per_vertex()
.unsigned int DoFTools::max_dofs_per_vertex | ( | const hp::DoFHandler< dim, spacedim > & | dh | ) |
Maximal number of degrees of freedom on a vertex.
This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.
dh.get_fe_collection().max_dofs_per_vertex()
.unsigned int DoFTools::n_components | ( | const DoFHandler< dim, spacedim > & | dh | ) |
Number of vector components in the finite element object used by this DoFHandler.
This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.
dh.get_fe_collection().n_components()
.unsigned int DoFTools::n_components | ( | const hp::DoFHandler< dim, spacedim > & | dh | ) |
Number of vector components in the finite element object used by this DoFHandler.
This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.
dh.get_fe_collection().n_components()
.bool DoFTools::fe_is_primitive | ( | const DoFHandler< dim, spacedim > & | dh | ) |
Find out whether the first FiniteElement used by this DoFHandler is primitive or not.
This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.
dh.get_fe(0).is_primitive()
.bool DoFTools::fe_is_primitive | ( | const hp::DoFHandler< dim, spacedim > & | dh | ) |
Find out whether the first FiniteElement used by this DoFHandler is primitive or not.
This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.
dh.get_fe(0).is_primitive()
.void DoFTools::make_periodicity_constraints | ( | const FaceIterator & | face_1, |
const typename identity< FaceIterator >::type & | face_2, | ||
::ConstraintMatrix & | constraint_matrix, | ||
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>() |
||
) |
Insert the (algebraic) constraints due to periodic boundary conditions into a ConstraintMatrix constraint_matrix
.
Given a pair of not necessarily active boundary faces face_1
and face_2
, this functions constrains all DoFs associated with the boundary described by face_1
to the respective DoFs of the boundary described by face_2
. More precisely:
If face_1
and face_2
are both active faces it adds the DoFs of face_1
to the list of constrained DoFs in constraint_matrix
and adds entries to constrain them to the corresponding values of the DoFs on face_2
. This happens on a purely algebraic level, meaning, the global DoF with (local face) index i
on face_1
gets constraint to the DoF with (local face) index i
on face_2
(possibly corrected for orientation, see below).
Otherwise, if face_1
and face_2
are not active faces, this function loops recursively over the children of face_1
and face_2
. If only one of the two faces is active, then we recursively iterate over the children of the non-active ones and make sure that the solution function on the refined side equals that on the non-refined face in much the same way as we enforce hanging node constraints at places where differently refined cells come together. (However, unlike hanging nodes, we do not enforce the requirement that there be only a difference of one refinement level between the two sides of the domain you would like to be periodic).
This routine only constrains DoFs that are not already constrained. If this routine encounters a DoF that already is constrained (for instance by Dirichlet boundary conditions), the old setting of the constraint (dofs the entry is constrained to, inhomogeneities) is kept and nothing happens.
The flags in the component_mask
(see GlossComponentMask) denote which components of the finite element space shall be constrained with periodic boundary conditions. If it is left as specified by the default value all components are constrained. If it is different from the default value, it is assumed that the number of entries equals the number of components of the finite element. This can be used to enforce periodicity in only one variable in a system of equations.
face_orientation
, face_flip
and face_rotation
describe an orientation that should be applied to face_1
prior to matching and constraining DoFs. This has nothing to do with the actual orientation of the given faces in their respective cells (which for boundary faces is always the default) but instead how you want to see periodicity to be enforced. For example, by using these flags, you can enforce a condition of the kind \(u(0,y)=u(1,1-y)\) (i.e., a Moebius band) or in 3d a twisted torus. More precisely, these flags match local face DoF indices in the following manner:
In 2d: face_orientation
must always be true
, face_rotation
is always false
, and face_flip has the meaning of line_flip
; this implies e.g. for Q1
:
And similarly for the case of Q1 in 3d:
Optionally a matrix matrix
along with a std::vector first_vector_components
can be specified that describes how DoFs on face_1
should be modified prior to constraining to the DoFs of face_2
. Here, two declarations are possible: If the std::vector first_vector_components
is non empty the matrix is interpreted as a dim
\(\times\) dim
rotation matrix that is applied to all vector valued blocks listed in first_vector_components
of the FESystem. If first_vector_components
is empty the matrix is interpreted as an interpolation matrix with size no_face_dofs \(\times\) no_face_dofs.
This function makes sure that identity constraints don't create cycles in constraint_matrix
.
Detailed information can be found in the see Glossary entry on periodic boundary conditions.
Definition at line 2199 of file dof_tools_constraints.cc.
void DoFTools::make_periodicity_constraints | ( | const std::vector< GridTools::PeriodicFacePair< typename DoFHandlerType::cell_iterator > > & | periodic_faces, |
::ConstraintMatrix & | constraint_matrix, | ||
const ComponentMask & | component_mask = ComponentMask() , |
||
const std::vector< unsigned int > & | first_vector_components = std::vector<unsigned int>() |
||
) |
Insert the (algebraic) constraints due to periodic boundary conditions into a ConstraintMatrix constraint_matrix
.
This is the main high level interface for above low level variant of make_periodicity_constraints(). It takes a std::vector periodic_faces
as argument and applies above make_periodicity_constraints() on each entry. periodic_faces
can be created by GridTools::collect_periodic_faces.
Definition at line 2424 of file dof_tools_constraints.cc.
void DoFTools::make_periodicity_constraints | ( | const DoFHandlerType & | dof_handler, |
const types::boundary_id | b_id1, | ||
const types::boundary_id | b_id2, | ||
const int | direction, | ||
::ConstraintMatrix & | constraint_matrix, | ||
const ComponentMask & | component_mask = ComponentMask() |
||
) |
Insert the (algebraic) constraints due to periodic boundary conditions into a ConstraintMatrix constraint_matrix
.
This function serves as a high level interface for the make_periodicity_constraints() function.
Define a 'first' boundary as all boundary faces having boundary_id b_id1
and a 'second' boundary consisting of all faces belonging to b_id2
.
This function tries to match all faces belonging to the first boundary with faces belonging to the second boundary with the help of orthogonal_equality().
If this matching is successful it constrains all DoFs associated with the 'first' boundary to the respective DoFs of the 'second' boundary respecting the relative orientation of the two faces.
Definition at line 2469 of file dof_tools_constraints.cc.
void DoFTools::make_periodicity_constraints | ( | const DoFHandlerType & | dof_handler, |
const types::boundary_id | b_id, | ||
const int | direction, | ||
::ConstraintMatrix & | constraint_matrix, | ||
const ComponentMask & | component_mask = ComponentMask() |
||
) |
This compatibility version of make_periodicity_constraints only works on grids with cells in standard orientation.
Instead of defining a 'first' and 'second' boundary with the help of two boundary_ids this function defines a 'left' boundary as all faces with local face index 2*dimension
and boundary indicator b_id
and, similarly, a 'right' boundary consisting of all face with local face index 2*dimension+1
and boundary indicator b_id
.
Definition at line 2500 of file dof_tools_constraints.cc.
void DoFTools::extract_hanging_node_dofs | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
std::vector< bool > & | selected_dofs | ||
) |
Select all dofs that will be constrained by interface constraints, i.e. all hanging nodes.
The size of selected_dofs
shall equal dof_handler.n_dofs()
. Previous contents of this array or overwritten.
In case of a parallel::shared::Triangulation or a parallel::distributed::Triangulation only locally relevant dofs are considered. Note that the vector returned through the second argument still has size dof_handler.n_dofs()
. Consequently, it can be very large for large parallel computations – in fact, it may be too large to store on each processor. In that case, you may want to choose the variant of this function that returns an IndexSet object.
Definition at line 981 of file dof_tools.cc.
IndexSet DoFTools::extract_hanging_node_dofs | ( | const DoFHandler< dim, spacedim > & | dof_handler | ) |
Same as above but return the selected DoFs as IndexSet. In particular, for parallel::Triangulation objects this function should be preferred.
Definition at line 997 of file dof_tools.cc.
void DoFTools::extract_dofs | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
const ComponentMask & | component_mask, | ||
std::vector< bool > & | selected_dofs | ||
) |
Extract the indices of the degrees of freedom belonging to certain vector components of a vector-valued finite element. The component_mask
defines which components or blocks of an FESystem are to be extracted from the DoFHandler dof
. The entries in the output array selected_dofs
corresponding to degrees of freedom belonging to these components are then flagged true
, while all others are set to false
.
If the finite element under consideration is not primitive, i.e., some or all of its shape functions are non-zero in more than one vector component (which holds, for example, for FE_Nedelec or FE_RaviartThomas elements), then shape functions cannot be associated with a single vector component. In this case, if one shape vector component of this element is flagged in component_mask
(see GlossComponentMask), then this is equivalent to selecting all vector components corresponding to this non-primitive base element.
[in] | dof_handler | The DoFHandler whose enumerated degrees of freedom are to be filtered by this function. |
[in] | component_mask | A mask that states which components you want to select. The size of this mask must be compatible with the number of components in the FiniteElement used by the dof_handler . See the glossary entry on component masks for more information. |
[out] | selected_dofs | A vector that will hold true or false values for each degree of freedom depending on whether or not it corresponds to a vector component selected by the mask above. The size of this array must equal DoFHandler::n_locally_owned_dofs(), which for sequential computations of course equals DoFHandler::n_dofs(). The previous contents of this array are overwritten. |
Definition at line 390 of file dof_tools.cc.
void DoFTools::extract_dofs | ( | const hp::DoFHandler< dim, spacedim > & | dof_handler, |
const ComponentMask & | component_mask, | ||
std::vector< bool > & | selected_dofs | ||
) |
The same function as above, but for a hp::DoFHandler.
Definition at line 436 of file dof_tools.cc.
void DoFTools::extract_dofs | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
const BlockMask & | block_mask, | ||
std::vector< bool > & | selected_dofs | ||
) |
This function is the equivalent to the DoFTools::extract_dofs() functions above except that the selection of which degrees of freedom to extract is not done based on components (see GlossComponent) but instead based on whether they are part of a particular block (see GlossBlock). Consequently, the second argument is not a ComponentMask but a BlockMask object.
[in] | dof_handler | The DoFHandler whose enumerated degrees of freedom are to be filtered by this function. |
[in] | block_mask | A mask that states which blocks you want to select. The size of this mask must be compatible with the number of blocks in the FiniteElement used by the dof_handler . See the glossary entry on block masks for more information. |
[out] | selected_dofs | A vector that will hold true or false values for each degree of freedom depending on whether or not it corresponds to a vector block selected by the mask above. The size of this array must equal DoFHandler::n_locally_owned_dofs(), which for sequential computations of course equals DoFHandler::n_dofs(). The previous contents of this array are overwritten. |
Definition at line 481 of file dof_tools.cc.
void DoFTools::extract_dofs | ( | const hp::DoFHandler< dim, spacedim > & | dof_handler, |
const BlockMask & | block_mask, | ||
std::vector< bool > & | selected_dofs | ||
) |
The same function as above, but for a hp::DoFHandler.
Definition at line 494 of file dof_tools.cc.
void DoFTools::extract_level_dofs | ( | const unsigned int | level, |
const DoFHandlerType & | dof, | ||
const ComponentMask & | component_mask, | ||
std::vector< bool > & | selected_dofs | ||
) |
Do the same thing as the corresponding extract_dofs() function for one level of a multi-grid DoF numbering.
Definition at line 507 of file dof_tools.cc.
void DoFTools::extract_level_dofs | ( | const unsigned int | level, |
const DoFHandlerType & | dof, | ||
const BlockMask & | component_mask, | ||
std::vector< bool > & | selected_dofs | ||
) |
Do the same thing as the corresponding extract_dofs() function for one level of a multi-grid DoF numbering.
Definition at line 559 of file dof_tools.cc.
void DoFTools::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>() |
||
) |
Extract all degrees of freedom which are at the boundary and belong to specified components of the solution. The function returns its results in the last non-default-valued parameter which contains true
if a degree of freedom is at the boundary and belongs to one of the selected components, and false
otherwise. The function is used in step-15.
By specifying the boundary_id
variable, you can select which boundary indicators the faces have to have on which the degrees of freedom are located that shall be extracted. If it is an empty list, then all boundary indicators are accepted.
The size of component_mask
(see GlossComponentMask) shall equal the number of components in the finite element used by dof
. The size of selected_dofs
shall equal dof_handler.n_dofs()
. Previous contents of this array or overwritten.
Using the usual convention, if a shape function is non-zero in more than one component (i.e. it is non-primitive), then the element in the component mask is used that corresponds to the first non-zero components. Elements in the mask corresponding to later components are ignored.
selected_dofs
has to have a length equal to all global degrees of freedom. Consequently, this does not scale to very large problems. If you need the functionality of this function for parallel triangulations, then you need to use the other DoFTools::extract_boundary_dofs function.[in] | dof_handler | The object that describes which degrees of freedom live on which cell |
[in] | component_mask | A mask denoting the vector components of the finite element that should be considered (see also GlossComponentMask). |
[out] | selected_dofs | A vector of booleans that is returned and for which an element will be true if the corresponding index is a degree of freedom that is located on the boundary (and correspond to the selected vector components and boundary indicators, depending on the values of the component_mask and boundary_ids arguments). |
[in] | boundary_ids | If empty, this function extracts the indices of the degrees of freedom for all parts of the boundary. If it is a non- empty list, then the function only considers boundary faces with the boundary indicators listed in this argument. |
Definition at line 573 of file dof_tools.cc.
void DoFTools::extract_boundary_dofs | ( | const DoFHandlerType & | dof_handler, |
const ComponentMask & | component_mask, | ||
IndexSet & | selected_dofs, | ||
const std::set< types::boundary_id > & | boundary_ids = std::set<types::boundary_id>() |
||
) |
This function does the same as the previous one but it returns its result as an IndexSet rather than a std::vector<bool>. Thus, it can also be called for DoFHandler objects that are defined on parallel::distributed::Triangulation objects.
selected_dofs
index set will contain only those degrees of freedom on the boundary that belong to the locally relevant set (see locally relevant DoFs).[in] | dof_handler | The object that describes which degrees of freedom live on which cell |
[in] | component_mask | A mask denoting the vector components of the finite element that should be considered (see also GlossComponentMask). |
[out] | selected_dofs | The IndexSet object that is returned and that will contain the indices of degrees of freedom that are located on the boundary (and correspond to the selected vector components and boundary indicators, depending on the values of the component_mask and boundary_ids arguments). |
[in] | boundary_ids | If empty, this function extracts the indices of the degrees of freedom for all parts of the boundary. If it is a non- empty list, then the function only considers boundary faces with the boundary indicators listed in this argument. |
Definition at line 599 of file dof_tools.cc.
void DoFTools::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>() |
||
) |
This function is similar to the extract_boundary_dofs() function but it extracts those degrees of freedom whose shape functions are nonzero on at least part of the selected boundary. For continuous elements, this is exactly the set of shape functions whose degrees of freedom are defined on boundary faces. On the other hand, if the finite element in used is a discontinuous element, all degrees of freedom are defined in the inside of cells and consequently none would be boundary degrees of freedom. Several of those would have shape functions that are nonzero on the boundary, however. This function therefore extracts all those for which the FiniteElement::has_support_on_face function says that it is nonzero on any face on one of the selected boundary parts.
Definition at line 707 of file dof_tools.cc.
IndexSet DoFTools::extract_dofs_with_support_contained_within | ( | const DoFHandlerType & | dof_handler, |
const std::function< bool(const typename DoFHandlerType::active_cell_iterator &)> & | predicate, | ||
const ConstraintMatrix & | cm = ConstraintMatrix() |
||
) |
Extract all indices of shape functions such that their support is entirely contained within the cells for which the predicate
is true
. The result is returned as an IndexSet.
Consider the following FE space where predicate returns true
for all cells on the left half of the domain:
This functions will return the union of all DoF indices on those cells minus DoF 11, 13, 2 and 0; the result will be [9,10], 12, [14,38]
. In the image above the returned DoFs are separated from the rest by the red line
Essentially, the question this functions answers is the following: Given a subdomain with associated DoFs, what is the largest subset of these DoFs that are allowed to be non-zero such that after calling ConstraintMatrix::distribute() the resulting solution vector will have support only within the given domain. Here cm
is the ConstraintMatrix containing hanging nodes constraints.
In case of parallel::distributed::Triangulation predicate
will be called only for locally owned and ghost cells. The resulting index set may contain DoFs that are associated with the locally owned or ghost cells, but are not owned by the current MPI core.
Definition at line 788 of file dof_tools.cc.
void DoFTools::extract_constant_modes | ( | const DoFHandlerType & | dof_handler, |
const ComponentMask & | component_mask, | ||
std::vector< std::vector< bool > > & | constant_modes | ||
) |
Extract a vector that represents the constant modes of the DoFHandler for the components chosen by component_mask
(see GlossComponentMask). The constant modes on a discretization are the null space of a Laplace operator on the selected components with Neumann boundary conditions applied. The null space is a necessary ingredient for obtaining a good AMG preconditioner when using the class TrilinosWrappers::PreconditionAMG. Since the ML AMG package only works on algebraic properties of the respective matrix, it has no chance to detect whether the matrix comes from a scalar or a vector valued problem. However, a near null space supplies exactly the needed information about the components placement of vector components within the matrix. The null space (or rather, the constant modes) is provided by the finite element underlying the given DoFHandler and for most elements, the null space will consist of as many vectors as there are true arguments in component_mask
(see GlossComponentMask), each of which will be one in one vector component and zero in all others. However, the representation of the constant function for e.g. FE_DGP is different (the first component on each element one, all other components zero), and some scalar elements may even have two constant modes (FE_Q_DG0). Therefore, we store this object in a vector of vectors, where the outer vector contains the collection of the actual constant modes on the DoFHandler. Each inner vector has as many components as there are (locally owned) degrees of freedom in the selected components. Note that any matrix associated with this null space must have been constructed using the same component_mask
argument, since the numbering of DoFs is done relative to the selected dofs, not to all dofs.
The main reason for this program is the use of the null space with the AMG preconditioner.
Definition at line 1175 of file dof_tools.cc.
void DoFTools::extract_subdomain_dofs | ( | const DoFHandlerType & | dof_handler, |
const types::subdomain_id | subdomain_id, | ||
std::vector< bool > & | selected_dofs | ||
) |
Flag all those degrees of freedom which are on cells with the given subdomain id. Note that DoFs on faces can belong to cells with differing subdomain ids, so the sets of flagged degrees of freedom are not mutually exclusive for different subdomain ids.
If you want to get a unique association of degree of freedom with subdomains, use the get_subdomain_association
function.
Definition at line 1006 of file dof_tools.cc.
void DoFTools::extract_locally_owned_dofs | ( | const DoFHandlerType & | dof_handler, |
IndexSet & | dof_set | ||
) |
Extract the set of global DoF indices that are owned by the current processor. For regular DoFHandler objects, this set is the complete set with all DoF indices. In either case, it equals what DoFHandler::locally_owned_dofs() returns.
Definition at line 1039 of file dof_tools.cc.
void DoFTools::extract_locally_active_dofs | ( | const DoFHandlerType & | dof_handler, |
IndexSet & | dof_set | ||
) |
Extract the set of global DoF indices that are active on the current DoFHandler. For regular DoFHandlers, these are all DoF indices, but for DoFHandler objects built on parallel::distributed::Triangulation this set is a superset of DoFHandler::locally_owned_dofs() and contains all DoF indices that live on all locally owned cells (including on the interface to ghost cells). However, it does not contain the DoF indices that are exclusively defined on ghost or artificial cells (see the glossary).
The degrees of freedom identified by this function equal those obtained from the dof_indices_with_subdomain_association() function when called with the locally owned subdomain id.
Definition at line 1051 of file dof_tools.cc.
void DoFTools::extract_locally_relevant_dofs | ( | const DoFHandlerType & | dof_handler, |
IndexSet & | dof_set | ||
) |
Extract the set of global DoF indices that are active on the current DoFHandler. For regular DoFHandlers, these are all DoF indices, but for DoFHandler objects built on parallel::distributed::Triangulation this set is the union of DoFHandler::locally_owned_dofs() and the DoF indices on all ghost cells. In essence, it is the DoF indices on all cells that are not artificial (see the glossary).
Definition at line 1087 of file dof_tools.cc.
std::vector< IndexSet > DoFTools::locally_owned_dofs_per_subdomain | ( | const DoFHandlerType & | dof_handler | ) |
For each processor, determine the set of locally owned degrees of freedom as an IndexSet. This function then returns a vector of index sets, where the vector has size equal to the number of MPI processes that participate in the DoF handler object.
The function can be used for objects of type Triangulation or parallel::shared::Triangulation. It will not work for objects of type parallel::distributed::Triangulation since for such triangulations we do not have information about all cells of the triangulation available locally, and consequently can not say anything definitive about the degrees of freedom active on other processors' locally owned cells.
Definition at line 1276 of file dof_tools.cc.
std::vector< IndexSet > DoFTools::locally_relevant_dofs_per_subdomain | ( | const DoFHandlerType & | dof_handler | ) |
For each processor, determine the set of locally relevant degrees of freedom as an IndexSet. This function then returns a vector of index sets, where the vector has size equal to the number of MPI processes that participate in the DoF handler object.
The function can be used for objects of type Triangulation or parallel::shared::Triangulation. It will not work for objects of type parallel::distributed::Triangulation since for such triangulations we do not have information about all cells of the triangulation available locally, and consequently can not say anything definitive about the degrees of freedom active on other processors' locally owned cells.
Definition at line 1367 of file dof_tools.cc.
void DoFTools::extract_locally_relevant_level_dofs | ( | const DoFHandlerType & | dof_handler, |
const unsigned int | level, | ||
IndexSet & | dof_set | ||
) |
Same as extract_locally_relevant_dofs() but for multigrid DoFs for the given level
.
Definition at line 1127 of file dof_tools.cc.
void DoFTools::get_subdomain_association | ( | const DoFHandlerType & | dof_handler, |
std::vector< types::subdomain_id > & | subdomain | ||
) |
For each degree of freedom, return in the output array to which subdomain (as given by the cell->subdomain_id()
function) it belongs. The output array is supposed to have the right size already when calling this function.
Note that degrees of freedom associated with faces, edges, and vertices may be associated with multiple subdomains if they are sitting on partition boundaries. In these cases, we assign them to the process with the smaller subdomain id. This may lead to different numbers of degrees of freedom in partitions, even if the number of cells is perfectly equidistributed. While this is regrettable, it is not a problem in practice since the number of degrees of freedom on partition boundaries is asymptotically vanishing as we refine the mesh as long as the number of partitions is kept constant.
This function returns the association of each DoF with one subdomain. If you are looking for the association of each cell with a subdomain, either query the cell->subdomain_id()
function, or use the GridTools::get_subdomain_association
function.
Note that this function is of questionable use for DoFHandler objects built on parallel::distributed::Triangulation since in that case ownership of individual degrees of freedom by MPI processes is controlled by the DoF handler object, not based on some geometric algorithm in conjunction with subdomain id. In particular, the degrees of freedom identified by the functions in this namespace as associated with a subdomain are not the same the DoFHandler class identifies as those it owns.
Definition at line 1430 of file dof_tools.cc.
unsigned int DoFTools::count_dofs_with_subdomain_association | ( | const DoFHandlerType & | dof_handler, |
const types::subdomain_id | subdomain | ||
) |
Count how many degrees of freedom are uniquely associated with the given subdomain
index.
Note that there may be rare cases where cells with the given subdomain
index exist, but none of its degrees of freedom are actually associated with it. In that case, the returned value will be zero.
This function will generate an exception if there are no cells with the given subdomain
index.
This function returns the number of DoFs associated with one subdomain. If you are looking for the association of cells with this subdomain, use the GridTools::count_cells_with_subdomain_association
function.
Note that this function is of questionable use for DoFHandler objects built on parallel::distributed::Triangulation since in that case ownership of individual degrees of freedom by MPI processes is controlled by the DoF handler object, not based on some geometric algorithm in conjunction with subdomain id. In particular, the degrees of freedom identified by the functions in this namespace as associated with a subdomain are not the same the DoFHandler class identifies as those it owns.
Definition at line 1517 of file dof_tools.cc.
void DoFTools::count_dofs_with_subdomain_association | ( | const DoFHandlerType & | dof_handler, |
const types::subdomain_id | subdomain, | ||
std::vector< unsigned int > & | n_dofs_on_subdomain | ||
) |
Count how many degrees of freedom are uniquely associated with the given subdomain
index.
This function does what the previous one does except that it splits the result among the vector components of the finite element in use by the DoFHandler object. The last argument (which must have a length equal to the number of vector components) will therefore store how many degrees of freedom of each vector component are associated with the given subdomain.
Note that this function is of questionable use for DoFHandler objects built on parallel::distributed::Triangulation since in that case ownership of individual degrees of freedom by MPI processes is controlled by the DoF handler object, not based on some geometric algorithm in conjunction with subdomain id. In particular, the degrees of freedom identified by the functions in this namespace as associated with a subdomain are not the same the DoFHandler class identifies as those it owns.
Definition at line 1590 of file dof_tools.cc.
IndexSet DoFTools::dof_indices_with_subdomain_association | ( | const DoFHandlerType & | dof_handler, |
const types::subdomain_id | subdomain | ||
) |
Return a set of indices that denotes the degrees of freedom that live on the given subdomain, i.e. that are on cells owned by the current processor. Note that this includes the ones that this subdomain "owns" (i.e. the ones for which get_subdomain_association() returns a value equal to the subdomain given here and that are selected by the DoFHandler::locally_owned_dofs() function) but also all of those that sit on the boundary between the given subdomain and other subdomain. In essence, degrees of freedom that sit on boundaries between subdomain will be in the index sets returned by this function for more than one subdomain.
Note that this function is of questionable use for DoFHandler objects built on parallel::distributed::Triangulation since in that case ownership of individual degrees of freedom by MPI processes is controlled by the DoF handler object, not based on some geometric algorithm in conjunction with subdomain id. In particular, the degrees of freedom identified by the functions in this namespace as associated with a subdomain are not the same the DoFHandler class identifies as those it owns.
Definition at line 1532 of file dof_tools.cc.
std::vector< types::global_dof_index > DoFTools::get_dofs_on_patch | ( | const std::vector< typename DoFHandlerType::active_cell_iterator > & | patch | ) |
Return the set of degrees of freedom that live on a set of cells (i.e., a patch) described by the argument.
Patches are often used in defining error estimators that require the solution of a local problem on the patch surrounding each of the cells of the mesh. You can get a list of cells that form the patch around a given cell using GridTools::get_patch_around_cell(). While DoFTools::count_dofs_on_patch() can be used to determine the size of these local problems, so that one can assemble the local system and then solve it, it is still necessary to provide a mapping between the global indices of the degrees of freedom that live on the patch and a local enumeration. This function provides such a local enumeration by returning the set of degrees of freedom that live on the patch.
Since this set is returned in the form of a std::vector, one can also think of it as a mapping
where i
is an index into the returned vector (i.e., a the local index of a degree of freedom on the patch) and global_dof_index
is the global index of a degree of freedom located on the patch. The array returned has size equal to DoFTools::count_dofs_on_patch().
DoFHandlerType | A type that is either DoFHandler or hp::DoFHandler. In C++, the compiler can not determine the type of DoFHandlerType from the function call. You need to specify it as an explicit template argument following the function name. |
patch | A collection of cells within an object of type DoFHandlerType |
Definition at line 2594 of file dof_tools.cc.
void DoFTools::make_cell_patches | ( | SparsityPattern & | block_list, |
const DoFHandler< dim, spacedim > & | dof_handler, | ||
const unsigned int | level, | ||
const std::vector< bool > & | selected_dofs = std::vector<bool>() , |
||
const types::global_dof_index | offset = 0 |
||
) |
Creates a sparsity pattern, which lists the degrees of freedom associated to each cell on the given level. This pattern can be used in RelaxationBlock classes as block list for additive and multiplicative Schwarz methods.
The row index in this pattern is the cell index resulting from standard iteration through a level of the Triangulation. For a parallel::distributed::Triangulation, only locally owned cells are entered.
The sparsity pattern is resized in this function to contain as many rows as there are locally owned cells on a given level, as many columns as there are degrees of freedom on this level.
selected_dofs
is a vector indexed by the local degrees of freedom on a cell. If it is used, only such dofs are entered into the block list which are selected. This allows for instance the exclusion of components or of dofs on the boundary.
Definition at line 2256 of file dof_tools.cc.
std::vector< unsigned int > DoFTools::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 |
||
) |
Create an incidence matrix that for every vertex on a given level of a multilevel DoFHandler flags which degrees of freedom are associated with the adjacent cells. This data structure is a matrix with as many rows as there are vertices on a given level, as many columns as there are degrees of freedom on this level, and entries that are either true or false. This data structure is conveniently represented by a SparsityPattern object. The sparsity pattern may be empty when entering this function and will be reinitialized to the correct size.
The function has some boolean arguments (listed below) controlling details of the generated patches. The default settings are those for Arnold-Falk-Winther type smoothers for divergence and curl conforming finite elements with essential boundary conditions. Other applications are possible, in particular changing boundary_patches
for non- essential boundary conditions.
This function returns the vertex_mapping
, that contains the mapping from the vertex indices to the block indices of the block_list
. For vertices that do not lead to a vertex patch, the entry in vertex_mapping
contains the value invalid_unsigned_int
. If invert_vertex_mapping
is set to true
, then the vertex_mapping
is inverted such that it contains the mapping from the block indices to the corresponding vertex indices.
block_list
: the SparsityPattern into which the patches will be stored.dof_handler
: the multilevel dof handler providing the topology operated on.interior_dofs_only
: for each patch of cells around a vertex, collect only the interior degrees of freedom of the patch and disregard those on the boundary of the patch. This is for instance the setting for smoothers of Arnold-Falk-Winther type.boundary_patches
: include patches around vertices at the boundary of the domain. If not, only patches around interior vertices will be generated.level_boundary_patches
: same for refinement edges towards coarser cells.single_cell_patches
: if not true, patches containing a single cell are eliminated.invert_vertex_mapping
: if true, then the return value contains one vertex index for each block; if false, then the return value contains one block index or invalid_unsigned_int
for each vertex. Definition at line 2405 of file dof_tools.cc.
std::vector< unsigned int > DoFTools::make_vertex_patches | ( | SparsityPattern & | block_list, |
const DoFHandlerType & | dof_handler, | ||
const unsigned int | level, | ||
const BlockMask & | exclude_boundary_dofs = BlockMask() , |
||
const bool | boundary_patches = false , |
||
const bool | level_boundary_patches = false , |
||
const bool | single_cell_patches = false , |
||
const bool | invert_vertex_mapping = false |
||
) |
Same as above but allows boundary dofs on blocks to be excluded individually.
This is helpful if you want to use, for example, Taylor Hood elements as it allows you to not include the boundary DoFs for the velocity block on the patches while also letting you include the boundary DoFs for the pressure block.
For each patch of cells around a vertex, collect all of the interior degrees of freedom of the patch and disregard those on the boundary of the patch if the boolean value for the corresponding block in the BlockMask of exclude_boundary_dofs
is false.
Definition at line 2428 of file dof_tools.cc.
void DoFTools::make_child_patches | ( | SparsityPattern & | block_list, |
const DoFHandlerType & | dof_handler, | ||
const unsigned int | level, | ||
const bool | interior_dofs_only, | ||
const bool | boundary_dofs = false |
||
) |
Create an incidence matrix that for every cell on a given level of a multilevel DoFHandler flags which degrees of freedom are associated with children of this cell. This data structure is conveniently represented by a SparsityPattern object.
The function thus creates a sparsity pattern which in each row (with rows corresponding to the cells on this level) lists the degrees of freedom associated to the cells that are the children of this cell. The DoF indices used here are level dof indices of a multilevel hierarchy, i.e., they may be associated with children that are not themselves active. The sparsity pattern may be empty when entering this function and will be reinitialized to the correct size.
The function has some boolean arguments (listed below) controlling details of the generated patches. The default settings are those for Arnold-Falk-Winther type smoothers for divergence and curl conforming finite elements with essential boundary conditions. Other applications are possible, in particular changing boundary_dofs
for non- essential boundary conditions.
block_list
: the SparsityPattern into which the patches will be stored.dof_handler
: The multilevel dof handler providing the topology operated on.interior_dofs_only
: for each patch of cells around a vertex, collect only the interior degrees of freedom of the patch and disregard those on the boundary of the patch. This is for instance the setting for smoothers of Arnold-Falk-Winther type.boundary_dofs
: include degrees of freedom, which would have excluded by interior_dofs_only
, but are lying on the boundary of the domain, and thus need smoothing. This parameter has no effect if interior_dofs_only
is false. Definition at line 2341 of file dof_tools.cc.
void DoFTools::make_single_patch | ( | SparsityPattern & | block_list, |
const DoFHandlerType & | dof_handler, | ||
const unsigned int | level, | ||
const bool | interior_dofs_only = false |
||
) |
Create a block list with only a single patch, which in turn contains all degrees of freedom on the given level.
This function is mostly a closure on level 0 for functions like make_child_patches() and make_vertex_patches(), which may produce an empty patch list.
block_list
: the SparsityPattern into which the patches will be stored.dof_handler
: The multilevel dof handler providing the topology operated on.level
The grid level used for building the list.interior_dofs_only
: if true, exclude degrees of freedom on the boundary of the domain. Definition at line 2298 of file dof_tools.cc.
void DoFTools::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 = std::vector<unsigned int>() |
||
) |
Count how many degrees of freedom out of the total number belong to each component. If the number of components the finite element has is one (i.e. you only have one scalar variable), then the number in this component obviously equals the total number of degrees of freedom. Otherwise, the sum of the DoFs in all the components needs to equal the total number.
However, the last statement does not hold true if the finite element is not primitive, i.e. some or all of its shape functions are non-zero in more than one vector component. This applies, for example, to the Nedelec or Raviart-Thomas elements. In this case, a degree of freedom is counted in each component in which it is non-zero, so that the sum mentioned above is greater than the total number of degrees of freedom.
This behavior can be switched off by the optional parameter vector_valued_once
. If this is true
, the number of components of a nonprimitive vector valued element is collected only in the first component. All other components will have a count of zero.
The additional optional argument target_component
allows for a re- sorting and grouping of components. To this end, it contains for each component the component number it shall be counted as. Having the same number entered several times sums up several components as the same. One of the applications of this argument is when you want to form block matrices and vectors, but want to pack several components into the same block (for example, when you have dim
velocities and one pressure, to put all velocities into one block, and the pressure into another).
The result is returned in dofs_per_component
. Note that the size of dofs_per_component
needs to be enough to hold all the indices specified in target_component
. If this is not the case, an assertion is thrown. The indices not targeted by target_components are left untouched.
Definition at line 1749 of file dof_tools.cc.
void DoFTools::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>() |
||
) |
Count the degrees of freedom in each block. This function is similar to count_dofs_per_component(), with the difference that the counting is done by blocks. See blocks in the glossary for details. Again the vectors are assumed to have the correct size before calling this function. If this is not the case, an assertion is thrown.
This function is used in the step-22, step-31, and step-32 tutorial programs.
Definition at line 1836 of file dof_tools.cc.
void DoFTools::get_active_fe_indices | ( | const DoFHandlerType & | dof_handler, |
std::vector< unsigned int > & | active_fe_indices | ||
) |
For each active cell of a DoFHandler or hp::DoFHandler, extract the active finite element index and fill the vector given as second argument. This vector is assumed to have as many entries as there are active cells.
For non-hp DoFHandler objects given as first argument, the returned vector will consist of only zeros, indicating that all cells use the same finite element. For a hp::DoFHandler, the values may be different, though.
Definition at line 1262 of file dof_tools.cc.
unsigned int DoFTools::count_dofs_on_patch | ( | const std::vector< typename DoFHandlerType::active_cell_iterator > & | patch | ) |
Count how many degrees of freedom live on a set of cells (i.e., a patch) described by the argument.
Patches are often used in defining error estimators that require the solution of a local problem on the patch surrounding each of the cells of the mesh. You can get a list of cells that form the patch around a given cell using GridTools::get_patch_around_cell(). This function is then useful in setting up the size of the linear system used to solve the local problem on the patch around a cell. The function DoFTools::get_dofs_on_patch() will then help to make the connection between global degrees of freedom and the local ones.
DoFHandlerType | A type that is either DoFHandler or hp::DoFHandler. In C++, the compiler can not determine the type of DoFHandlerType from the function call. You need to specify it as an explicit template argument following the function name. |
patch | A collection of cells within an object of type DoFHandlerType |
Definition at line 2566 of file dof_tools.cc.
void DoFTools::map_dof_to_boundary_indices | ( | const DoFHandlerType & | dof_handler, |
std::vector< types::global_dof_index > & | mapping | ||
) |
Create a mapping from degree of freedom indices to the index of that degree of freedom on the boundary. After this operation, mapping[dof]
gives the index of the degree of freedom with global number dof
in the list of degrees of freedom on the boundary. If the degree of freedom requested is not on the boundary, the value of mapping[dof]
is numbers::invalid_dof_index. This function is mainly used when setting up matrices and vectors on the boundary from the trial functions, which have global numbers, while the matrices and vectors use numbers of the trial functions local to the boundary.
Prior content of mapping
is deleted.
Definition at line 1916 of file dof_tools.cc.
void DoFTools::map_dof_to_boundary_indices | ( | const DoFHandlerType & | dof_handler, |
const std::set< types::boundary_id > & | boundary_ids, | ||
std::vector< types::global_dof_index > & | mapping | ||
) |
Same as the previous function, except that only those parts of the boundary are considered for which the boundary indicator is listed in the second argument.
See the general doc of this class for more information.
Definition at line 1955 of file dof_tools.cc.
void DoFTools::map_dofs_to_support_points | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandler< dim, spacedim > & | dof_handler, | ||
std::vector< Point< spacedim > > & | support_points | ||
) |
Return a list of support points (see this glossary entry) for all the degrees of freedom handled by this DoF handler object. This function, of course, only works if the finite element object used by the DoF handler object actually provides support points, i.e. no edge elements or the like. Otherwise, an exception is thrown.
Definition at line 2077 of file dof_tools.cc.
void DoFTools::map_dofs_to_support_points | ( | const ::hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::DoFHandler< dim, spacedim > & | dof_handler, | ||
std::vector< Point< spacedim > > & | support_points | ||
) |
Same as the previous function but for the hp case.
void DoFTools::map_dofs_to_support_points | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandler< dim, spacedim > & | dof_handler, | ||
std::map< types::global_dof_index, Point< spacedim > > & | support_points | ||
) |
This function is a version of the above map_dofs_to_support_points function that doesn't simply return a vector of support points (see this glossary entry) with one entry for each global degree of freedom, but instead a map that maps from the DoFs index to its location. The point of this function is that it is also usable in cases where the DoFHandler is based on a parallel::distributed::Triangulation object. In such cases, each processor will not be able to determine the support point location of all DoFs, and worse no processor may be able to hold a vector that would contain the locations of all DoFs even if they were known. As a consequence, this function constructs a map from those DoFs for which we can know the locations (namely, those DoFs that are locally relevant (see locally relevant DoFs) to their locations.
For non-distributed triangulations, the map returned as support_points
is of course dense, i.e., every DoF is to be found in it.
mapping | The mapping from the reference cell to the real cell on which DoFs are defined. |
dof_handler | The object that describes which DoF indices live on which cell of the triangulation. |
support_points | A map that for every locally relevant DoF index contains the corresponding location in real space coordinates. Previous content of this object is deleted in this function. |
Definition at line 2123 of file dof_tools.cc.
void DoFTools::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 | ||
) |
Same as the previous function but for the hp case.
void DoFTools::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 | ||
) |
This is the opposite function to the one above. It generates a map where the keys are the support points of the degrees of freedom, while the values are the DoF indices. For a definition of support points, see this glossary entry.
Since there is no natural order in the space of points (except for the 1d case), you have to provide a map with an explicitly specified comparator object. This function is therefore templatized on the comparator object. Previous content of the map object is deleted in this function.
Just as with the function above, it is assumed that the finite element in use here actually supports the notion of support points of all its components.
void DoFTools::distribute_cell_to_dof_vector | ( | const DoFHandlerType & | dof_handler, |
const Vector< Number > & | cell_data, | ||
Vector< double > & | dof_data, | ||
const unsigned int | component = 0 |
||
) |
Take a vector of values which live on cells (e.g. an error per cell) and distribute it to the dofs in such a way that a finite element field results, which can then be further processed, e.g. for output. You should note that the resulting field will not be continuous at hanging nodes. This can, however, easily be arranged by calling the appropriate distribute
function of a ConstraintMatrix object created for this DoFHandler object, after the vector has been fully assembled.
It is assumed that the number of elements in cell_data
equals the number of active cells and that the number of elements in dof_data
equals dof_handler.n_dofs()
.
Note that the input vector may be a vector of any data type as long as it is convertible to double
. The output vector, being a data vector on a DoF handler, always consists of elements of type double
.
In case the finite element used by this DoFHandler consists of more than one component, you need to specify which component in the output vector should be used to store the finite element field in; the default is zero (no other value is allowed if the finite element consists only of one component). All other components of the vector remain untouched, i.e. their contents are not changed.
This function cannot be used if the finite element in use has shape functions that are non-zero in more than one vector component (in deal.II speak: they are non-primitive).
Definition at line 308 of file dof_tools.cc.
void DoFTools::write_gnuplot_dof_support_point_info | ( | std::ostream & | out, |
const std::map< types::global_dof_index, Point< spacedim > > & | support_points | ||
) |
Generate text output readable by gnuplot with point data based on the given map support_points
. For each support point location, a string label containing a list of all DoFs from the map is generated. The map can be generated with a call to map_dofs_to_support_points() and is useful to visualize location and global numbering of unknowns.
An example for the format of each line in the output is:
where x, y, and z (present only in corresponding dimension) are the coordinates of the support point, followed by a list of DoF numbers.
The points with labels can be plotted as follows in gnuplot:
Examples (this also includes the grid written separately using GridOut):
To generate the mesh and the support point information in a single gnuplot file, use code similar to
and from within gnuplot execute the following command:
Alternatively, the following gnuplot script will generate a png file when executed as gnuplot gnuplot.gpl
on the command line:
Definition at line 2156 of file dof_tools.cc.