Reference documentation for deal.II version GIT relicensing-468-ge3ce94fd06 2024-04-23 15:40:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Namespaces | Functions
MGTools Namespace Reference

Namespaces

namespace  internal
 

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<int dim, int spacedim, typename number = double>
void make_sparsity_pattern (const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity, const unsigned int level, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true)
 
template<int dim, int spacedim, typename number = double>
void make_flux_sparsity_pattern (const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity, const unsigned int level, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true)
 
template<int dim, int spacedim>
void make_flux_sparsity_pattern_edge (const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity, const unsigned int level)
 
template<int dim, int spacedim>
void make_flux_sparsity_pattern (const DoFHandler< dim, spacedim > &dof, SparsityPatternBase &sparsity, const unsigned int level, const Table< 2, DoFTools::Coupling > &int_mask, const Table< 2, DoFTools::Coupling > &flux_mask)
 
template<int dim, int spacedim>
void make_flux_sparsity_pattern_edge (const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity, const unsigned int level, const Table< 2, DoFTools::Coupling > &flux_mask)
 
template<int dim, int spacedim>
void make_interface_sparsity_pattern (const DoFHandler< dim, spacedim > &dof_handler, const MGConstrainedDoFs &mg_constrained_dofs, SparsityPatternBase &sparsity, const unsigned int level)
 
template<int dim, int spacedim>
void count_dofs_per_block (const DoFHandler< dim, spacedim > &dof_handler, std::vector< std::vector< types::global_dof_index > > &dofs_per_block, std::vector< unsigned int > target_block={})
 
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={})
 
template<int dim, int spacedim>
void make_boundary_list (const DoFHandler< dim, spacedim > &mg_dof, const std::map< types::boundary_id, const Function< spacedim > * > &function_map, std::vector< std::set< types::global_dof_index > > &boundary_indices, const ComponentMask &component_mask={})
 
template<int dim, int spacedim>
void make_boundary_list (const DoFHandler< dim, spacedim > &mg_dof, const std::map< types::boundary_id, const Function< spacedim > * > &function_map, std::vector< IndexSet > &boundary_indices, const ComponentMask &component_mask={})
 
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={})
 
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)
 
template<int dim, int spacedim>
std::vector< types::global_dof_indexlocal_workload (const Triangulation< dim, spacedim > &tria)
 
template<int dim, int spacedim>
std::vector< types::global_dof_indexlocal_workload (const std::vector< std::shared_ptr< const Triangulation< dim, spacedim > > > &trias)
 
template<int dim, int spacedim>
double workload_imbalance (const Triangulation< dim, spacedim > &tria)
 
template<int dim, int spacedim>
double workload_imbalance (const std::vector< std::shared_ptr< const Triangulation< dim, spacedim > > > &trias)
 
template<int dim, int spacedim>
std::vector< std::pair< types::global_dof_index, types::global_dof_index > > local_vertical_communication_cost (const Triangulation< dim, spacedim > &tria)
 
template<int dim, int spacedim>
std::vector< std::pair< types::global_dof_index, types::global_dof_index > > local_vertical_communication_cost (const std::vector< std::shared_ptr< const Triangulation< dim, spacedim > > > &trias)
 
template<int dim, int spacedim>
double vertical_communication_efficiency (const Triangulation< dim, spacedim > &tria)
 
template<int dim, int spacedim>
double vertical_communication_efficiency (const std::vector< std::shared_ptr< const Triangulation< dim, spacedim > > > &trias)
 
template<>
void compute_row_length_vector (const DoFHandler< 1, 1 > &, const unsigned int, std::vector< unsigned int > &, const DoFTools::Coupling)
 
template<>
void compute_row_length_vector (const DoFHandler< 1, 1 > &, const unsigned int, std::vector< unsigned int > &, const Table< 2, DoFTools::Coupling > &, const Table< 2, DoFTools::Coupling > &)
 
template<>
void compute_row_length_vector (const DoFHandler< 1, 2 > &, const unsigned int, std::vector< unsigned int > &, const DoFTools::Coupling)
 
template<>
void compute_row_length_vector (const DoFHandler< 1, 2 > &, const unsigned int, std::vector< unsigned int > &, const Table< 2, DoFTools::Coupling > &, const Table< 2, DoFTools::Coupling > &)
 

Detailed Description

addtogroup mg 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.

Function Documentation

◆ compute_row_length_vector() [1/6]

template<int dim, int spacedim>
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 106 of file mg_tools.cc.

◆ compute_row_length_vector() [2/6]

template<int dim, int spacedim>
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 291 of file mg_tools.cc.

◆ make_sparsity_pattern()

template<int dim, int spacedim, typename number = double>
void MGTools::make_sparsity_pattern ( const DoFHandler< dim, spacedim > &  dof_handler,
SparsityPatternBase sparsity,
const unsigned int  level,
const AffineConstraints< number > &  constraints = {},
const bool  keep_constrained_dofs = true 
)

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, either using SparsityPattern::compress() or by copying to a SparsityPattern.

The optional AffineConstraints argument allows to define constraints of the level matrices like Dirichlet boundary conditions. Note that there is need to consider hanging nodes on the typical level matrices, since only one level is considered. See DoFTools::make_sparsity_pattern() for more details about the arguments.

Definition at line 575 of file mg_tools.cc.

◆ make_flux_sparsity_pattern() [1/2]

template<int dim, int spacedim, typename number = double>
void MGTools::make_flux_sparsity_pattern ( const DoFHandler< dim, spacedim > &  dof_handler,
SparsityPatternBase sparsity,
const unsigned int  level,
const AffineConstraints< number > &  constraints = {},
const bool  keep_constrained_dofs = true 
)

Make a sparsity pattern including fluxes of discontinuous Galerkin methods.

See also
make_sparsity_pattern and DoFTools

Definition at line 605 of file mg_tools.cc.

◆ make_flux_sparsity_pattern_edge() [1/2]

template<int dim, int spacedim>
void MGTools::make_flux_sparsity_pattern_edge ( const DoFHandler< dim, spacedim > &  dof_handler,
SparsityPatternBase 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.

make_flux_sparsity_pattern()

Definition at line 679 of file mg_tools.cc.

◆ make_flux_sparsity_pattern() [2/2]

template<int dim, int spacedim>
void MGTools::make_flux_sparsity_pattern ( const DoFHandler< dim, spacedim > &  dof,
SparsityPatternBase 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 741 of file mg_tools.cc.

◆ make_flux_sparsity_pattern_edge() [2/2]

template<int dim, int spacedim>
void MGTools::make_flux_sparsity_pattern_edge ( const DoFHandler< dim, spacedim > &  dof_handler,
SparsityPatternBase 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.

make_flux_sparsity_pattern()

Definition at line 919 of file mg_tools.cc.

◆ make_interface_sparsity_pattern()

template<int dim, int spacedim>
void MGTools::make_interface_sparsity_pattern ( const DoFHandler< dim, spacedim > &  dof_handler,
const MGConstrainedDoFs mg_constrained_dofs,
SparsityPatternBase 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 1013 of file mg_tools.cc.

◆ count_dofs_per_block()

template<int dim, int spacedim>
void MGTools::count_dofs_per_block ( const DoFHandler< dim, spacedim > &  dof_handler,
std::vector< std::vector< types::global_dof_index > > &  dofs_per_block,
std::vector< unsigned int target_block = {} 
)

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 1156 of file mg_tools.cc.

◆ count_dofs_per_component()

template<int dim, int spacedim>
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 = {} 
)

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 1054 of file mg_tools.cc.

◆ make_boundary_list() [1/3]

template<int dim, int spacedim>
void MGTools::make_boundary_list ( const DoFHandler< dim, spacedim > &  mg_dof,
const std::map< types::boundary_id, const Function< spacedim > * > &  function_map,
std::vector< std::set< types::global_dof_index > > &  boundary_indices,
const ComponentMask component_mask = {} 
)

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 1239 of file mg_tools.cc.

◆ make_boundary_list() [2/3]

template<int dim, int spacedim>
void MGTools::make_boundary_list ( const DoFHandler< dim, spacedim > &  mg_dof,
const std::map< types::boundary_id, const Function< spacedim > * > &  function_map,
std::vector< IndexSet > &  boundary_indices,
const ComponentMask component_mask = {} 
)

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 1264 of file mg_tools.cc.

◆ make_boundary_list() [3/3]

template<int dim, int spacedim>
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 = {} 
)

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 1285 of file mg_tools.cc.

◆ extract_inner_interface_dofs()

template<int dim, int spacedim>
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 1443 of file mg_tools.cc.

◆ max_level_for_coarse_mesh()

template<int dim, int spacedim>
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 1530 of file mg_tools.cc.

◆ local_workload() [1/2]

template<int dim, int spacedim>
std::vector< types::global_dof_index > MGTools::local_workload ( const Triangulation< dim, spacedim > &  tria)

This function returns the local work of each level, i.e., the number of locally-owned cells on each multigrid level.

Note
This function requires that parallel::TriangulationBase::is_multilevel_hierarchy_constructed() is true, which can be controlled by setting the construct_multigrid_hierarchy flag when constructing the Triangulation.

Definition at line 1628 of file mg_tools.cc.

◆ local_workload() [2/2]

template<int dim, int spacedim>
std::vector< types::global_dof_index > MGTools::local_workload ( const std::vector< std::shared_ptr< const Triangulation< dim, spacedim > > > &  trias)

Similar to the above function but for a vector of triangulations as created e.g. by MGTransferGlobalCoarseningTools::create_geometric_coarsening_sequence(). This returns the number of locally-owned cells on each of the triangulations.

Definition at line 1654 of file mg_tools.cc.

◆ workload_imbalance() [1/2]

template<int dim, int spacedim>
double MGTools::workload_imbalance ( const Triangulation< dim, spacedim > &  tria)

Return the imbalance of the parallel distribution of the multigrid mesh hierarchy, based on the local workload provided by local_workload(). Ideally this value is equal to 1 (every processor owns the same number of cells on each level, approximately true for most globally refined meshes). Values greater than 1 estimate the slowdown one should see in a geometric multigrid v-cycle as compared with the same computation on a perfectly distributed mesh hierarchy.

Note
This function is a collective MPI call between all ranks of the Triangulation and therefore needs to be called from all ranks.
This function requires that parallel::TriangulationBase::is_multilevel_hierarchy_constructed() is true, which can be controlled by setting the construct_multigrid_hierarchy flag when constructing the Triangulation.

Definition at line 1674 of file mg_tools.cc.

◆ workload_imbalance() [2/2]

template<int dim, int spacedim>
double MGTools::workload_imbalance ( const std::vector< std::shared_ptr< const Triangulation< dim, spacedim > > > &  trias)

Similar to the above function but for a vector of triangulations as created e.g. by MGTransferGlobalCoarseningTools::create_geometric_coarsening_sequence().

Definition at line 1684 of file mg_tools.cc.

◆ local_vertical_communication_cost() [1/2]

template<int dim, int spacedim>
std::vector< std::pair< types::global_dof_index, types::global_dof_index > > MGTools::local_vertical_communication_cost ( const Triangulation< dim, spacedim > &  tria)

Return the vertical communication cost between levels. The returned vector contains for each level the number of cells that have the same owning process as their corresponding coarse cell (parent) and the number of cells that have not the same owning process as their corresponding coarse cell.

Definition at line 1696 of file mg_tools.cc.

◆ local_vertical_communication_cost() [2/2]

template<int dim, int spacedim>
std::vector< std::pair< types::global_dof_index, types::global_dof_index > > MGTools::local_vertical_communication_cost ( const std::vector< std::shared_ptr< const Triangulation< dim, spacedim > > > &  trias)

Similar to the above function but for a vector of triangulations as created e.g. by MGTransferGlobalCoarseningTools::create_geometric_coarsening_sequence().

Note
This function is a collective MPI call between all ranks of the Triangulation and therefore needs to be called from all ranks.

Definition at line 1730 of file mg_tools.cc.

◆ vertical_communication_efficiency() [1/2]

template<int dim, int spacedim>
double MGTools::vertical_communication_efficiency ( const Triangulation< dim, spacedim > &  tria)

Share of fine cells that have the same owning process as their corresponding coarse cell (parent). This quantity gives an indication on the efficiency of a multigrid transfer operator and on how much data has to be sent around. A small number indicates that most of the data has to be completely permuted, involving a large volume of communication.

Note
This function is a collective MPI call between all ranks of the Triangulation and therefore needs to be called from all ranks.

Definition at line 1805 of file mg_tools.cc.

◆ vertical_communication_efficiency() [2/2]

template<int dim, int spacedim>
double MGTools::vertical_communication_efficiency ( const std::vector< std::shared_ptr< const Triangulation< dim, spacedim > > > &  trias)

Similar to the above function but for a vector of triangulations as created e.g. by MGTransferGlobalCoarseningTools::create_geometric_coarsening_sequence().

Definition at line 1815 of file mg_tools.cc.

◆ compute_row_length_vector() [3/6]

template<>
void MGTools::compute_row_length_vector ( const DoFHandler< 1, 1 > &  ,
const unsigned int  ,
std::vector< unsigned int > &  ,
const DoFTools::Coupling   
)

Definition at line 56 of file mg_tools.cc.

◆ compute_row_length_vector() [4/6]

template<>
void MGTools::compute_row_length_vector ( const DoFHandler< 1, 1 > &  ,
const unsigned int  ,
std::vector< unsigned int > &  ,
const Table< 2, DoFTools::Coupling > &  ,
const Table< 2, DoFTools::Coupling > &   
)

Definition at line 68 of file mg_tools.cc.

◆ compute_row_length_vector() [5/6]

template<>
void MGTools::compute_row_length_vector ( const DoFHandler< 1, 2 > &  ,
const unsigned int  ,
std::vector< unsigned int > &  ,
const DoFTools::Coupling   
)

Definition at line 81 of file mg_tools.cc.

◆ compute_row_length_vector() [6/6]

template<>
void MGTools::compute_row_length_vector ( const DoFHandler< 1, 2 > &  ,
const unsigned int  ,
std::vector< unsigned int > &  ,
const Table< 2, DoFTools::Coupling > &  ,
const Table< 2, DoFTools::Coupling > &   
)

Definition at line 92 of file mg_tools.cc.