FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
FEAT::Solver::AmaVanka< Matrix_, Filter_ > Class Template Reference

Additive Macro-wise Matrix-based Vanka preconditioner/smoother. More...

#include <amavanka.hpp>

Inheritance diagram for FEAT::Solver::AmaVanka< Matrix_, Filter_ >:
FEAT::Solver::SolverBase< Matrix_::VectorTypeL > FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >

Public Types

typedef Solver::SolverBase< typename Matrix_::VectorTypeL > BaseClass
 our base-class More...
 
typedef Matrix_::DataType DataType
 our data type More...
 
typedef Matrix_::IndexType IndexType
 our index type More...
 
typedef Matrix_::VectorTypeL VectorType
 our vector type More...
 

Public Member Functions

 AmaVanka (const Matrix_ &matrix, const Filter_ &filter, const DataType omega=DataType(1), const Index num_steps=Index(1))
 Constructor. More...
 
virtual Status apply (VectorType &vec_x, const VectorType &vec_b) override
 applies the preconditioner More...
 
std::size_t bytes () const
 Returns the total number of bytes currently allocated in this object. More...
 
void clear_macro_dofs ()
 Clears the macro dofs graphs. More...
 
bool compare (const AmaVanka *other) const
 
std::size_t data_size () const
 Returns the total data size used by the AmaVanka smoother. More...
 
virtual void done ()
 Finalization method. More...
 
virtual void done_numeric ()
 Numeric finalization method. More...
 
virtual void done_symbolic () override
 Releases the symbolic factorization data. More...
 
virtual void init ()
 Initialization method. More...
 
virtual void init_numeric () override
 Performs numeric factorization. More...
 
virtual void init_symbolic () override
 Performs symbolic factorization. More...
 
virtual String name () const override
 Returns the name of the solver. More...
 
void push_macro_dofs (Adjacency::Graph &&dofs)
 Pushes the dofs-at-macro graph of the next block. More...
 
void reset_timings ()
 Resets the internal stop watches for time measurement. More...
 
void set_num_steps (Index num_steps)
 Sets the number of smoothing steps. More...
 
void set_omega (DataType omega)
 Sets the damping parameter omega. More...
 
void set_skip_singular (bool skip_sing)
 Sets whether singular macros are to be skipped. More...
 
double time_apply () const
 Returns the total accumulated time for the solver application. More...
 
double time_init_numeric () const
 Returns the total accumulated time for numeric initialization. More...
 
double time_init_symbolic () const
 Returns the total accumulated time for symbolic initialization. More...
 

Protected Types

typedef Intern::AmaVankaMatrixHelper< Matrix_ >::VankaMatrix VankaMatrixType
 the type of our Vanka matrix More...
 

Protected Attributes

bool _auto_macros
 deduct macro dofs automatically? More...
 
std::vector< Adjacency::Graph_dof_macros
 
const Filter_ & _filter
 the system filter More...
 
std::vector< Adjacency::Graph_macro_dofs
 the DOF-macro graphs More...
 
std::vector< int > _macro_mask
 the macro mask More...
 
const Matrix_ & _matrix
 the system matrix More...
 
Index _num_steps
 number of steps More...
 
int _num_threads
 number of threads for numeric factorization More...
 
DataType _omega
 damping parameter More...
 
bool _skip_singular
 skip singular macros? More...
 
VankaMatrixType _vanka
 the Vanka preconditioner matrix More...
 
VectorType _vec_c
 temporary vectors More...
 
VectorType _vec_d
 
StopWatch watch_apply
 
StopWatch watch_init_numeric
 
StopWatch watch_init_symbolic
 

Detailed Description

template<typename Matrix_, typename Filter_>
class FEAT::Solver::AmaVanka< Matrix_, Filter_ >

Additive Macro-wise Matrix-based Vanka preconditioner/smoother.

This class implements an additive macro-wise Vanka smoother, which stores its pre-computed operator as a sparse matrix, so that each application of the Vanka smoother consists of only one sparse matrix-vector multiplication.

This class supports a whole batch of different matrix types:

Attention
In contrast to the Solver::Vanka class, one must specify the macros for the Vanka smoother for each finite element space which is used to discretize the individual blocks of the solution space by using the push_macro_dofs() function. The only exception is for the LAFEM::SaddlePointMatrix, as for this matrix type, there exists an automatic macro deduction algorithm, which will compute the element-wise macro graphs based on the sparsity patterns of the B and D matrices, resp.

Example: Assume that your solution space is the usual velocity-pressure space pair used in a Stokes system, then you need to add the dof-mapping graphs of both the velocity and the pressure space to the smoother in the following way:

auto vanka = Solver::new_amavanka(system_matrix, system_filter);
vanka->push_macro_dofs(Space::DofMappingRenderer::render(space_velocity));
vanka->push_macro_dofs(Space::DofMappingRenderer::render(space_pressure));
static Adjacency::Graph render(const Space_ &space)
Renders the dof-mapping of a space into an adjacency graph.
std::shared_ptr< AmaVanka< Matrix_, Filter_ > > new_amavanka(const Matrix_ &matrix, const Filter_ &filter, typename Matrix_::DataType omega=typename Matrix_::DataType(1), Index num_steps=Index(1))
Creates a new AmaVanka smoother object.
Definition: amavanka.hpp:497

If your solution space has other components (such as e.g. temperature, density or stress), then you have to add the rendered dof-mappings of those spaces in the correct order as well.

Note
In the case of a LAFEM::SaddlePointMatrix, the smoother implemented in this class is mathematically equivalent to a Solver::Vanka smoother of type Solver::VankaType::block_full_add.
Author
Peter Zajac

Definition at line 63 of file amavanka.hpp.

Member Typedef Documentation

◆ BaseClass

template<typename Matrix_ , typename Filter_ >
typedef Solver::SolverBase<typename Matrix_::VectorTypeL> FEAT::Solver::AmaVanka< Matrix_, Filter_ >::BaseClass

our base-class

Definition at line 68 of file amavanka.hpp.

◆ DataType

template<typename Matrix_ , typename Filter_ >
typedef Matrix_::DataType FEAT::Solver::AmaVanka< Matrix_, Filter_ >::DataType

our data type

Definition at line 71 of file amavanka.hpp.

◆ IndexType

template<typename Matrix_ , typename Filter_ >
typedef Matrix_::IndexType FEAT::Solver::AmaVanka< Matrix_, Filter_ >::IndexType

our index type

Definition at line 73 of file amavanka.hpp.

◆ VankaMatrixType

template<typename Matrix_ , typename Filter_ >
typedef Intern::AmaVankaMatrixHelper<Matrix_>::VankaMatrix FEAT::Solver::AmaVanka< Matrix_, Filter_ >::VankaMatrixType
protected

the type of our Vanka matrix

Definition at line 79 of file amavanka.hpp.

◆ VectorType

template<typename Matrix_ , typename Filter_ >
typedef Matrix_::VectorTypeL FEAT::Solver::AmaVanka< Matrix_, Filter_ >::VectorType

our vector type

Definition at line 75 of file amavanka.hpp.

Constructor & Destructor Documentation

◆ AmaVanka()

template<typename Matrix_ , typename Filter_ >
FEAT::Solver::AmaVanka< Matrix_, Filter_ >::AmaVanka ( const Matrix_ &  matrix,
const Filter_ &  filter,
const DataType  omega = DataType(1),
const Index  num_steps = Index(1) 
)
inlineexplicit

Constructor.

Parameters
[in]matrixThe saddle-point system matrix.
[in]filterThe system filter.
[in]omegaThe damping parameter.
[in]num_stepsThe number of smoothing steps to be performed.

Definition at line 128 of file amavanka.hpp.

Member Function Documentation

◆ apply()

template<typename Matrix_ , typename Filter_ >
virtual Status FEAT::Solver::AmaVanka< Matrix_, Filter_ >::apply ( VectorType vec_x,
const VectorType vec_b 
)
inlineoverridevirtual

◆ bytes()

template<typename Matrix_ , typename Filter_ >
std::size_t FEAT::Solver::AmaVanka< Matrix_, Filter_ >::bytes ( ) const
inline

Returns the total number of bytes currently allocated in this object.

Definition at line 212 of file amavanka.hpp.

References FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_macro_dofs, FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_vanka, and FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_vec_c.

◆ clear_macro_dofs()

template<typename Matrix_ , typename Filter_ >
void FEAT::Solver::AmaVanka< Matrix_, Filter_ >::clear_macro_dofs ( )
inline

Clears the macro dofs graphs.

Definition at line 143 of file amavanka.hpp.

References FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_auto_macros.

◆ compare()

template<typename Matrix_ , typename Filter_ >
bool FEAT::Solver::AmaVanka< Matrix_, Filter_ >::compare ( const AmaVanka< Matrix_, Filter_ > *  other) const
inline

Definition at line 472 of file amavanka.hpp.

◆ data_size()

template<typename Matrix_ , typename Filter_ >
std::size_t FEAT::Solver::AmaVanka< Matrix_, Filter_ >::data_size ( ) const
inline

Returns the total data size used by the AmaVanka smoother.

Returns
The total data size, i.e. the total number of floating point values used in the factorization.

Definition at line 230 of file amavanka.hpp.

References FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_vanka.

◆ done()

virtual void FEAT::Solver::SolverBase< Matrix_::VectorTypeL >::done ( )
inlinevirtualinherited

Finalization method.

This function performs both the symbolic and numeric finalization, i.e. it simply performs

this->done_numeric();
this->done_symbolic();

Definition at line 283 of file base.hpp.

◆ done_numeric()

virtual void FEAT::Solver::SolverBase< Matrix_::VectorTypeL >::done_numeric ( )
inlinevirtualinherited

Numeric finalization method.

This method is called to release any data allocated in the numeric initialization step.

Reimplemented in FEAT::Solver::GenericUmfpack< Matrix_ >.

Definition at line 246 of file base.hpp.

◆ done_symbolic()

template<typename Matrix_ , typename Filter_ >
virtual void FEAT::Solver::AmaVanka< Matrix_, Filter_ >::done_symbolic ( )
inlineoverridevirtual

◆ init()

virtual void FEAT::Solver::SolverBase< Matrix_::VectorTypeL >::init ( )
inlinevirtualinherited

Initialization method.

This function performs both the symbolic and numeric initialization, i.e. it simply performs

this->init_symbolic();
this->init_numeric();

Definition at line 268 of file base.hpp.

◆ init_numeric()

template<typename Matrix_ , typename Filter_ >
virtual void FEAT::Solver::AmaVanka< Matrix_, Filter_ >::init_numeric ( )
inlineoverridevirtual

◆ init_symbolic()

◆ name()

template<typename Matrix_ , typename Filter_ >
virtual String FEAT::Solver::AmaVanka< Matrix_, Filter_ >::name ( ) const
inlineoverridevirtual

Returns the name of the solver.

Implements FEAT::Solver::SolverBase< Matrix_::VectorTypeL >.

Reimplemented in FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >.

Definition at line 270 of file amavanka.hpp.

◆ push_macro_dofs()

template<typename Matrix_ , typename Filter_ >
void FEAT::Solver::AmaVanka< Matrix_, Filter_ >::push_macro_dofs ( Adjacency::Graph &&  dofs)
inline

Pushes the dofs-at-macro graph of the next block.

Parameters
[in]dofsThe dofs-at-macro graph of the block.

Definition at line 155 of file amavanka.hpp.

References FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_auto_macros, FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_macro_dofs, and XABORTM.

◆ reset_timings()

template<typename Matrix_ , typename Filter_ >
void FEAT::Solver::AmaVanka< Matrix_, Filter_ >::reset_timings ( )
inline

Resets the internal stop watches for time measurement.

Definition at line 238 of file amavanka.hpp.

References FEAT::StopWatch::reset().

◆ set_num_steps()

template<typename Matrix_ , typename Filter_ >
void FEAT::Solver::AmaVanka< Matrix_, Filter_ >::set_num_steps ( Index  num_steps)
inline

Sets the number of smoothing steps.

Parameters
[in]num_stepsThe number of smoothing steps to be performed.

Definition at line 180 of file amavanka.hpp.

References XASSERT.

◆ set_omega()

template<typename Matrix_ , typename Filter_ >
void FEAT::Solver::AmaVanka< Matrix_, Filter_ >::set_omega ( DataType  omega)
inline

Sets the damping parameter omega.

Parameters
[in]omegaThe damping parameter.

Definition at line 192 of file amavanka.hpp.

References XASSERT.

◆ set_skip_singular()

template<typename Matrix_ , typename Filter_ >
void FEAT::Solver::AmaVanka< Matrix_, Filter_ >::set_skip_singular ( bool  skip_sing)
inline

Sets whether singular macros are to be skipped.

Parameters
[in]skip_singSpecifies whether singular macros are to be skipped.

Definition at line 204 of file amavanka.hpp.

◆ time_apply()

template<typename Matrix_ , typename Filter_ >
double FEAT::Solver::AmaVanka< Matrix_, Filter_ >::time_apply ( ) const
inline

Returns the total accumulated time for the solver application.

Definition at line 264 of file amavanka.hpp.

References FEAT::StopWatch::elapsed().

◆ time_init_numeric()

template<typename Matrix_ , typename Filter_ >
double FEAT::Solver::AmaVanka< Matrix_, Filter_ >::time_init_numeric ( ) const
inline

Returns the total accumulated time for numeric initialization.

Definition at line 256 of file amavanka.hpp.

References FEAT::StopWatch::elapsed().

◆ time_init_symbolic()

template<typename Matrix_ , typename Filter_ >
double FEAT::Solver::AmaVanka< Matrix_, Filter_ >::time_init_symbolic ( ) const
inline

Returns the total accumulated time for symbolic initialization.

Definition at line 248 of file amavanka.hpp.

References FEAT::StopWatch::elapsed().

Member Data Documentation

◆ _auto_macros

template<typename Matrix_ , typename Filter_ >
bool FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_auto_macros
protected

◆ _dof_macros

template<typename Matrix_ , typename Filter_ >
std::vector<Adjacency::Graph> FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_dof_macros
protected

Definition at line 92 of file amavanka.hpp.

◆ _filter

template<typename Matrix_ , typename Filter_ >
const Filter_& FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_filter
protected

the system filter

Definition at line 84 of file amavanka.hpp.

◆ _macro_dofs

◆ _macro_mask

template<typename Matrix_ , typename Filter_ >
std::vector<int> FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_macro_mask
protected

◆ _matrix

template<typename Matrix_ , typename Filter_ >
const Matrix_& FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_matrix
protected

◆ _num_steps

template<typename Matrix_ , typename Filter_ >
Index FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_num_steps
protected

number of steps

Definition at line 96 of file amavanka.hpp.

Referenced by FEAT::Solver::AmaVanka< Matrix_, Filter_ >::apply().

◆ _num_threads

template<typename Matrix_ , typename Filter_ >
int FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_num_threads
protected

number of threads for numeric factorization

Definition at line 110 of file amavanka.hpp.

◆ _omega

◆ _skip_singular

◆ _vanka

◆ _vec_c

template<typename Matrix_ , typename Filter_ >
VectorType FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_vec_c
protected

temporary vectors

Definition at line 100 of file amavanka.hpp.

Referenced by FEAT::Solver::AmaVanka< Matrix_, Filter_ >::bytes().

◆ _vec_d

template<typename Matrix_ , typename Filter_ >
VectorType FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_vec_d
protected

Definition at line 100 of file amavanka.hpp.

◆ watch_apply

template<typename Matrix_ , typename Filter_ >
StopWatch FEAT::Solver::AmaVanka< Matrix_, Filter_ >::watch_apply
protected

Definition at line 107 of file amavanka.hpp.

◆ watch_init_numeric

template<typename Matrix_ , typename Filter_ >
StopWatch FEAT::Solver::AmaVanka< Matrix_, Filter_ >::watch_init_numeric
protected

Definition at line 105 of file amavanka.hpp.

◆ watch_init_symbolic

template<typename Matrix_ , typename Filter_ >
StopWatch FEAT::Solver::AmaVanka< Matrix_, Filter_ >::watch_init_symbolic
protected

Definition at line 103 of file amavanka.hpp.


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