Reference documentation for deal.II version 9.0.0
|
Functions | |
template<int dim, int spacedim> | |
void | compute_row_length_vector (const DoFHandler< dim, spacedim > &dofs, const unsigned int level, std::vector< unsigned int > &row_lengths, const DoFTools::Coupling flux_couplings=DoFTools::none) |
template<int dim, int spacedim> | |
void | compute_row_length_vector (const DoFHandler< dim, spacedim > &dofs, const unsigned int level, std::vector< unsigned int > &row_lengths, const Table< 2, DoFTools::Coupling > &couplings, const Table< 2, DoFTools::Coupling > &flux_couplings) |
template<typename DoFHandlerType , typename SparsityPatternType > | |
void | make_sparsity_pattern (const DoFHandlerType &dof_handler, SparsityPatternType &sparsity, const unsigned int level) |
template<int dim, typename SparsityPatternType , int spacedim> | |
void | make_flux_sparsity_pattern (const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity, const unsigned int level) |
template<int dim, typename SparsityPatternType , int spacedim> | |
void | make_flux_sparsity_pattern_edge (const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity, const unsigned int level) |
template<int dim, typename SparsityPatternType , int spacedim> | |
void | make_flux_sparsity_pattern (const DoFHandler< dim, spacedim > &dof, SparsityPatternType &sparsity, const unsigned int level, const Table< 2, DoFTools::Coupling > &int_mask, const Table< 2, DoFTools::Coupling > &flux_mask) |
template<int dim, typename SparsityPatternType , int spacedim> | |
void | make_flux_sparsity_pattern_edge (const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity, const unsigned int level, const Table< 2, DoFTools::Coupling > &flux_mask) |
template<typename DoFHandlerType , typename SparsityPatternType > | |
void | make_interface_sparsity_pattern (const DoFHandlerType &dof_handler, const MGConstrainedDoFs &mg_constrained_dofs, SparsityPatternType &sparsity, const unsigned int level) |
template<typename DoFHandlerType > | |
void | count_dofs_per_block (const DoFHandlerType &dof_handler, std::vector< std::vector< types::global_dof_index > > &dofs_per_block, std::vector< unsigned int > target_block=std::vector< unsigned int >()) |
template<int dim, int spacedim> | |
void | count_dofs_per_component (const DoFHandler< dim, spacedim > &mg_dof, std::vector< std::vector< types::global_dof_index > > &result, const bool only_once=false, std::vector< unsigned int > target_component=std::vector< unsigned int >()) |
template<int dim, int spacedim> | |
void | make_boundary_list (const DoFHandler< dim, spacedim > &mg_dof, const typename FunctionMap< dim >::type &function_map, std::vector< std::set< types::global_dof_index > > &boundary_indices, const ComponentMask &component_mask=ComponentMask()) |
template<int dim, int spacedim> | |
void | make_boundary_list (const DoFHandler< dim, spacedim > &mg_dof, const typename FunctionMap< dim >::type &function_map, std::vector< IndexSet > &boundary_indices, const ComponentMask &component_mask=ComponentMask()) |
template<int dim, int spacedim> | |
void | make_boundary_list (const DoFHandler< dim, spacedim > &mg_dof, const std::set< types::boundary_id > &boundary_ids, std::vector< IndexSet > &boundary_indices, const ComponentMask &component_mask=ComponentMask()) |
template<int dim, int spacedim> | |
void | extract_inner_interface_dofs (const DoFHandler< dim, spacedim > &mg_dof_handler, std::vector< IndexSet > &interface_dofs) |
template<int dim, int spacedim> | |
unsigned int | max_level_for_coarse_mesh (const Triangulation< dim, spacedim > &tria) |
This is a collection of functions operating on, and manipulating the numbers of degrees of freedom in a multilevel triangulation. It is similar in purpose and function to the DoFTools
namespace, but operates on levels of DoFHandler objects. See there and the documentation of the member functions for more information.
void MGTools::compute_row_length_vector | ( | const DoFHandler< dim, spacedim > & | dofs, |
const unsigned int | level, | ||
std::vector< unsigned int > & | row_lengths, | ||
const DoFTools::Coupling | flux_couplings = DoFTools::none |
||
) |
Compute row length vector for multilevel methods.
Definition at line 107 of file mg_tools.cc.
void MGTools::compute_row_length_vector | ( | const DoFHandler< dim, spacedim > & | dofs, |
const unsigned int | level, | ||
std::vector< unsigned int > & | row_lengths, | ||
const Table< 2, DoFTools::Coupling > & | couplings, | ||
const Table< 2, DoFTools::Coupling > & | flux_couplings | ||
) |
Compute row length vector for multilevel methods with optimization for block couplings.
Definition at line 280 of file mg_tools.cc.
void MGTools::make_sparsity_pattern | ( | const DoFHandlerType & | dof_handler, |
SparsityPatternType & | sparsity, | ||
const unsigned int | level | ||
) |
Write the sparsity structure of the matrix belonging to the specified level
. The sparsity pattern is not compressed, so before creating the actual matrix you have to compress the matrix yourself, using SparseMatrixStruct::compress()
.
There is no need to consider hanging nodes here, since only one level is considered.
Definition at line 540 of file mg_tools.cc.
void MGTools::make_flux_sparsity_pattern | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
SparsityPatternType & | sparsity, | ||
const unsigned int | level | ||
) |
Make a sparsity pattern including fluxes of discontinuous Galerkin methods.
Definition at line 573 of file mg_tools.cc.
void MGTools::make_flux_sparsity_pattern_edge | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
SparsityPatternType & | sparsity, | ||
const unsigned int | level | ||
) |
Create sparsity pattern for the fluxes at refinement edges. The matrix maps a function of the fine level space level
to the coarser space.
Definition at line 648 of file mg_tools.cc.
void MGTools::make_flux_sparsity_pattern | ( | const DoFHandler< dim, spacedim > & | dof, |
SparsityPatternType & | sparsity, | ||
const unsigned int | level, | ||
const Table< 2, DoFTools::Coupling > & | int_mask, | ||
const Table< 2, DoFTools::Coupling > & | flux_mask | ||
) |
This function does the same as the other with the same name, but it gets two additional coefficient matrices. A matrix entry will only be generated for two basis functions, if there is a non-zero entry linking their associated components in the coefficient matrix.
There is one matrix for couplings in a cell and one for the couplings occurring in fluxes.
Definition at line 716 of file mg_tools.cc.
void MGTools::make_flux_sparsity_pattern_edge | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
SparsityPatternType & | sparsity, | ||
const unsigned int | level, | ||
const Table< 2, DoFTools::Coupling > & | flux_mask | ||
) |
Create sparsity pattern for the fluxes at refinement edges. The matrix maps a function of the fine level space level
to the coarser space. This is the version restricting the pattern to the elements actually needed.
Definition at line 898 of file mg_tools.cc.
void MGTools::make_interface_sparsity_pattern | ( | const DoFHandlerType & | dof_handler, |
const MGConstrainedDoFs & | mg_constrained_dofs, | ||
SparsityPatternType & | sparsity, | ||
const unsigned int | level | ||
) |
Create sparsity pattern for interface_in/out matrices used in a multigrid computation. These matrices contain an entry representing the coupling of degrees of freedom on a refinement edge to those not on the refinement edge of a certain level.
Definition at line 988 of file mg_tools.cc.
void MGTools::count_dofs_per_block | ( | const DoFHandlerType & | dof_handler, |
std::vector< std::vector< types::global_dof_index > > & | dofs_per_block, | ||
std::vector< unsigned int > | target_block = std::vector<unsigned int>() |
||
) |
Count the dofs block-wise on each level.
Result is a vector containing for each level a vector containing the number of dofs for each block (access is result[level][block]
).
Definition at line 1125 of file mg_tools.cc.
void MGTools::count_dofs_per_component | ( | const DoFHandler< dim, spacedim > & | mg_dof, |
std::vector< std::vector< types::global_dof_index > > & | result, | ||
const bool | only_once = false , |
||
std::vector< unsigned int > | target_component = std::vector<unsigned int>() |
||
) |
Count the dofs component-wise on each level.
Result is a vector containing for each level a vector containing the number of dofs for each component (access is result[level][component]
).
Definition at line 1021 of file mg_tools.cc.
void MGTools::make_boundary_list | ( | const DoFHandler< dim, spacedim > & | mg_dof, |
const typename FunctionMap< dim >::type & | function_map, | ||
std::vector< std::set< types::global_dof_index > > & | boundary_indices, | ||
const ComponentMask & | component_mask = ComponentMask() |
||
) |
Generate a list of those degrees of freedom at the boundary of the domain that should be eliminated from the matrix because they will be constrained by Dirichlet boundary conditions.
This is the multilevel equivalent of VectorTools::interpolate_boundary_values, but since the multilevel method does not have its own right hand side, the function values returned by the function object that is part of the function_map argument are ignored.
boundary_indices
is a vector which on return contains all indices of degrees of freedom for each level that are at the part of the boundary identified by the function_map argument. Its length has to match the number of levels in the dof handler object.Previous content in boundary_indices
is not overwritten, but added to.
Definition at line 1208 of file mg_tools.cc.
void MGTools::make_boundary_list | ( | const DoFHandler< dim, spacedim > & | mg_dof, |
const typename FunctionMap< dim >::type & | function_map, | ||
std::vector< IndexSet > & | boundary_indices, | ||
const ComponentMask & | component_mask = ComponentMask() |
||
) |
The same function as above, but return an IndexSet rather than a std::set<unsigned int> on each level.
Previous content in boundary_indices
is not overwritten, but added to.
Definition at line 1234 of file mg_tools.cc.
void MGTools::make_boundary_list | ( | const DoFHandler< dim, spacedim > & | mg_dof, |
const std::set< types::boundary_id > & | boundary_ids, | ||
std::vector< IndexSet > & | boundary_indices, | ||
const ComponentMask & | component_mask = ComponentMask() |
||
) |
The same function as above, but return an IndexSet rather than a std::set<unsigned int> on each level and use a std::set of boundary_ids as input.
Previous content in boundary_indices
is not overwritten, but added to.
Definition at line 1255 of file mg_tools.cc.
void MGTools::extract_inner_interface_dofs | ( | const DoFHandler< dim, spacedim > & | mg_dof_handler, |
std::vector< IndexSet > & | interface_dofs | ||
) |
For each level in a multigrid hierarchy, produce an IndexSet that indicates which of the degrees of freedom are along interfaces of this level to cells that only exist on coarser levels.
Definition at line 1466 of file mg_tools.cc.
unsigned int MGTools::max_level_for_coarse_mesh | ( | const Triangulation< dim, spacedim > & | tria | ) |
Return the highest possible level that can be used as the coarsest level in a Multigrid computation, that is, the highest level in the hierarchy whose mesh covers the entire domain. This corresponds to the minimum level of a cell on the active mesh. Since each processor only has a local view of the mesh, each processor must call this function. Note that this is a global minimum over the entire mesh and therefore each processor will return the same value.
Definition at line 1555 of file mg_tools.cc.