ExaDG
|
#include <grid.h>
Public Member Functions | |
MultigridMappings (std::shared_ptr< dealii::Mapping< dim > > mapping, std::shared_ptr< dealii::Mapping< dim > > mapping_coarse_levels) | |
MultigridMappings (std::shared_ptr< MappingDoFVector< dim, Number > > mapping_dof_vector, unsigned int const degree_coarse_mappings) | |
void | initialize_coarse_mappings (Grid< dim > const &grid, unsigned int const n_h_levels) |
dealii::Mapping< dim > const & | get_mapping (unsigned int const h_level, unsigned int const n_h_levels) const |
Public Attributes | |
std::function< void(std::shared_ptr< dealii::Triangulation< dim > const > const &fine_triangulation, std::vector< std::shared_ptr< dealii::Triangulation< dim > const > > const &coarse_triangulations)> | lambda_initialize_coarse_mappings |
This class handles mappings for use in multigrid.
The lambda function initialize_coarse_mappings() can be overwritten in order to realize a user-specific construction of the mapping on coarser grids.
|
inline |
Use this constructor if the fine-level mapping is a "normal" dealii::Mapping.
|
inline |
Use this constructor if the fine-level mapping is of type ExaDG::MappingDoFVector.
|
inline |
Returns the dealii::Mapping for a given h_level of n_h_levels.
|
inline |
Initializes the multigrid mappings on coarse h levels.
std::function<void( std::shared_ptr<dealii::Triangulation<dim> const> const & fine_triangulation, std::vector<std::shared_ptr<dealii::Triangulation<dim> const>> const & coarse_triangulations)> ExaDG::MultigridMappings< dim, Number >::lambda_initialize_coarse_mappings |
This function creates coarse mappings of type MappingDoFVector for use in multigrid with h-transfer given a fine level mapping of type MappingDoFVector.
The default implementation uses an interpolation of the fine-level mapping to coarser grids. You can overwrite this function in order to realize a user-specific construction of the mapping on coarser grids.