ExaDG
Loading...
Searching...
No Matches
ExaDG::MultigridPreconditionerBase< dim, Number, MultigridNumber_ > Class Template Referenceabstract
Inheritance diagram for ExaDG::MultigridPreconditionerBase< dim, Number, MultigridNumber_ >:
ExaDG::PreconditionerBase< Number > ExaDG::Poisson::MultigridPreconditioner< dim, Number, 1 > ExaDG::ConvDiff::MultigridPreconditioner< dim, Number > ExaDG::IncNS::MultigridPreconditioner< dim, Number > ExaDG::IncNS::MultigridPreconditionerProjection< dim, Number > ExaDG::Poisson::MultigridPreconditioner< dim, Number, n_components > ExaDG::Structure::MultigridPreconditioner< dim, Number >

Public Types

typedef MultigridNumber_ MultigridNumber
Public Types inherited from ExaDG::PreconditionerBase< Number >
typedef dealii::LinearAlgebra::distributed::Vector< Number > VectorType

Public Member Functions

 MultigridPreconditionerBase (MPI_Comm const &comm)
void initialize (MultigridData const &data, std::shared_ptr< Grid< dim > const > grid, std::shared_ptr< MultigridMappings< dim, Number > > const multigrid_mappings, dealii::FiniteElement< dim > const &fe, bool const operator_is_singular, Map_DBC const &dirichlet_bc, Map_DBC_ComponentMask const &dirichlet_bc_component_mask, bool const initialize_preconditioners)
void update () override
void vmult (VectorType &dst, VectorType const &src) const override
unsigned int solve (VectorType &dst, VectorType const &src) const
virtual void apply_smoother_on_fine_level (VectorTypeMG &dst, VectorTypeMG const &src) const
std::shared_ptr< TimerTreeget_timings () const override
Public Member Functions inherited from ExaDG::PreconditionerBase< Number >
bool needs_update () const

Protected Types

typedef std::map< dealii::types::boundary_id, std::shared_ptr< dealii::Function< dim > > > Map_DBC
typedef std::map< dealii::types::boundary_id, dealii::ComponentMask > Map_DBC_ComponentMask
typedef std::vector< dealii::GridTools::PeriodicFacePair< typename dealii::Triangulation< dim >::cell_iterator > > PeriodicFacePairs
typedef dealii::LinearAlgebra::distributed::Vector< Number > VectorType
typedef dealii::LinearAlgebra::distributed::Vector< MultigridNumber > VectorTypeMG

Protected Member Functions

void initialize_mapping ()
virtual void initialize_matrix_free_objects ()
void update_matrix_free_objects ()
void update_smoothers ()
void update_coarse_solver ()
virtual void initialize_dof_handler_and_constraints (bool is_singular, unsigned int const n_components, Map_DBC const &dirichlet_bc, Map_DBC_ComponentMask const &dirichlet_bc_component_mask)
void do_initialize_dof_handler_and_constraints (bool is_singular, unsigned int const n_components, Map_DBC const &dirichlet_bc, Map_DBC_ComponentMask const &dirichlet_bc_component_mask, dealii::MGLevelObject< std::shared_ptr< dealii::DoFHandler< dim > const > > &dofhandlers, dealii::MGLevelObject< std::shared_ptr< dealii::AffineConstraints< MultigridNumber > > > &constraints)
virtual void initialize_transfer_operators ()
void do_initialize_transfer_operators (std::shared_ptr< MultigridTransfer< dim, MultigridNumber, VectorTypeMG > > &transfers, unsigned int const dof_index)
unsigned int get_number_of_levels () const
dealii::Mapping< dim > const & get_mapping (unsigned int const h_level) const
void for_all_levels (std::function< void(unsigned int const)> const &function_on_level)
void for_all_smoothing_levels (std::function< void(unsigned int const)> const &function_on_level)
void transfer_from_fine_to_coarse_levels (std::function< void(unsigned int const, unsigned int const)> const &levelwise_transfer)

Protected Attributes

std::shared_ptr< Grid< dim > const > grid
std::shared_ptr< MultigridMappings< dim, Number > > multigrid_mappings
dealii::MGLevelObject< std::shared_ptr< dealii::DoFHandler< dim > const > > dof_handlers
dealii::MGLevelObject< std::shared_ptr< dealii::AffineConstraints< MultigridNumber > > > constraints
dealii::MGLevelObject< std::shared_ptr< MatrixFreeData< dim, MultigridNumber > > > matrix_free_data_objects
dealii::MGLevelObject< std::shared_ptr< dealii::MatrixFree< dim, MultigridNumber > > > matrix_free_objects
dealii::MGLevelObject< std::shared_ptr< Operator > > operators
std::shared_ptr< MultigridTransfer< dim, MultigridNumber, VectorTypeMG > > transfers
std::vector< MGLevelInfolevel_info
Protected Attributes inherited from ExaDG::PreconditionerBase< Number >
bool update_needed

Member Function Documentation

◆ for_all_levels()

template<int dim, typename Number, typename MultigridNumber_ = float>
void ExaDG::MultigridPreconditionerBase< dim, Number, MultigridNumber_ >::for_all_levels ( std::function< void(unsigned int const)> const & function_on_level)
inlineprotected

This is a generic function allowing to loop over all multigrid levels (including the coarsest level). The operation to be performed on each level is passed as a lambda with argument level.

◆ for_all_smoothing_levels()

template<int dim, typename Number, typename MultigridNumber_ = float>
void ExaDG::MultigridPreconditionerBase< dim, Number, MultigridNumber_ >::for_all_smoothing_levels ( std::function< void(unsigned int const)> const & function_on_level)
inlineprotected

This is a generic function allowing to loop over all smoothing levels (excluding the coarsest level). The operation to be performed on each level is passed as a lambda with argument level.

◆ get_number_of_levels()

template<int dim, typename Number, typename MultigridNumber>
unsigned int ExaDG::MultigridPreconditionerBase< dim, Number, MultigridNumber >::get_number_of_levels ( ) const
protected

Returns the number of levels.

The number of levels includes the coarse level and the finer smoothing levels, i.e. n_levels = 1 if the multigrid preconditioner is a coarse-grid solve on the coarse level only.

◆ get_timings()

template<int dim, typename Number, typename MultigridNumber>
std::shared_ptr< TimerTree > ExaDG::MultigridPreconditionerBase< dim, Number, MultigridNumber >::get_timings ( ) const
overridevirtual

◆ transfer_from_fine_to_coarse_levels()

template<int dim, typename Number, typename MultigridNumber_ = float>
void ExaDG::MultigridPreconditionerBase< dim, Number, MultigridNumber_ >::transfer_from_fine_to_coarse_levels ( std::function< void(unsigned int const, unsigned int const)> const & levelwise_transfer)
inlineprotected

This is a generic function allowing to successively transfer information from the fine level to all coarser multigrid levels. The operation to be performed for a transfer between two successive levels is passed as a lambda with fine_level as the first argument and coarse_level as the second argument.

◆ update()

template<int dim, typename Number, typename MultigridNumber>
void ExaDG::MultigridPreconditionerBase< dim, Number, MultigridNumber >::update ( )
overridevirtual

◆ update_coarse_solver()

template<int dim, typename Number, typename MultigridNumber>
void ExaDG::MultigridPreconditionerBase< dim, Number, MultigridNumber >::update_coarse_solver ( )
protected

This function updates the coarse-grid solver. The prerequisite to call this function is that the coarse-grid operator has been updated.

◆ update_smoothers()

template<int dim, typename Number, typename MultigridNumber>
void ExaDG::MultigridPreconditionerBase< dim, Number, MultigridNumber >::update_smoothers ( )
protected

This function updates the smoother for all smoothing levels. The prerequisite to call this function is that the multigrid operators have been updated.

◆ vmult()

template<int dim, typename Number, typename MultigridNumber>
void ExaDG::MultigridPreconditionerBase< dim, Number, MultigridNumber >::vmult ( VectorType & dst,
VectorType const & src ) const
overridevirtual

The documentation for this class was generated from the following files: