|
FEAT 3
Finite Element Analysis Toolbox
|
Additive Macro-wise Matrix-based Vanka preconditioner/smoother. More...
#include <amavanka.hpp>
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 |
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:
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:
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.
Definition at line 63 of file amavanka.hpp.
| typedef Solver::SolverBase<typename Matrix_::VectorTypeL> FEAT::Solver::AmaVanka< Matrix_, Filter_ >::BaseClass |
our base-class
Definition at line 68 of file amavanka.hpp.
| typedef Matrix_::DataType FEAT::Solver::AmaVanka< Matrix_, Filter_ >::DataType |
our data type
Definition at line 71 of file amavanka.hpp.
| typedef Matrix_::IndexType FEAT::Solver::AmaVanka< Matrix_, Filter_ >::IndexType |
our index type
Definition at line 73 of file amavanka.hpp.
|
protected |
the type of our Vanka matrix
Definition at line 79 of file amavanka.hpp.
| typedef Matrix_::VectorTypeL FEAT::Solver::AmaVanka< Matrix_, Filter_ >::VectorType |
our vector type
Definition at line 75 of file amavanka.hpp.
|
inlineexplicit |
Constructor.
| [in] | matrix | The saddle-point system matrix. |
| [in] | filter | The system filter. |
| [in] | omega | The damping parameter. |
| [in] | num_steps | The number of smoothing steps to be performed. |
Definition at line 128 of file amavanka.hpp.
|
inlineoverridevirtual |
applies the preconditioner
Implements FEAT::Solver::SolverBase< Matrix_::VectorTypeL >.
Definition at line 444 of file amavanka.hpp.
References FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_num_steps, FEAT::StopWatch::start(), FEAT::StopWatch::stop(), and FEAT::Solver::success.
|
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.
|
inline |
Clears the macro dofs graphs.
Definition at line 143 of file amavanka.hpp.
References FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_auto_macros.
|
inline |
Definition at line 472 of file amavanka.hpp.
|
inline |
Returns the total data size used by the AmaVanka smoother.
Definition at line 230 of file amavanka.hpp.
References FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_vanka.
|
inlinevirtualinherited |
|
inlinevirtualinherited |
Numeric finalization method.
This method is called to release any data allocated in the numeric initialization step.
Reimplemented in FEAT::Solver::GenericUmfpack< Matrix_ >.
|
inlineoverridevirtual |
Releases the symbolic factorization data.
Reimplemented from FEAT::Solver::SolverBase< Matrix_::VectorTypeL >.
Reimplemented in FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >.
Definition at line 322 of file amavanka.hpp.
References FEAT::Solver::SolverBase< typename Matrix_::VectorTypeL >::done_symbolic().
Referenced by FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::done_symbolic().
|
inlinevirtualinherited |
|
inlineoverridevirtual |
Performs numeric factorization.
Reimplemented from FEAT::Solver::SolverBase< Matrix_::VectorTypeL >.
Reimplemented in FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >.
Definition at line 338 of file amavanka.hpp.
References FEAT::Solver::SolverBase< typename Matrix_::VectorTypeL >::init_numeric(), FEAT::Math::invert_matrix(), FEAT::StopWatch::start(), FEAT::StopWatch::stop(), and XASSERTM.
|
inlineoverridevirtual |
Performs symbolic factorization.
Reimplemented from FEAT::Solver::SolverBase< Matrix_::VectorTypeL >.
Reimplemented in FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >.
Definition at line 276 of file amavanka.hpp.
References FEAT::Solver::AmaVanka< Matrix_, Filter_ >::_macro_dofs, FEAT::Solver::SolverBase< typename Matrix_::VectorTypeL >::init_symbolic(), FEAT::StopWatch::start(), FEAT::StopWatch::stop(), FEAT::Adjacency::transpose, XABORTM, XASSERT, and XASSERTM.
Referenced by FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::init_symbolic().
|
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.
|
inline |
Pushes the dofs-at-macro graph of the next block.
| [in] | dofs | The 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.
|
inline |
Resets the internal stop watches for time measurement.
Definition at line 238 of file amavanka.hpp.
References FEAT::StopWatch::reset().
|
inline |
Sets the number of smoothing steps.
| [in] | num_steps | The number of smoothing steps to be performed. |
Definition at line 180 of file amavanka.hpp.
References XASSERT.
|
inline |
Sets the damping parameter omega.
| [in] | omega | The damping parameter. |
Definition at line 192 of file amavanka.hpp.
References XASSERT.
|
inline |
Sets whether singular macros are to be skipped.
| [in] | skip_sing | Specifies whether singular macros are to be skipped. |
Definition at line 204 of file amavanka.hpp.
|
inline |
Returns the total accumulated time for the solver application.
Definition at line 264 of file amavanka.hpp.
References FEAT::StopWatch::elapsed().
|
inline |
Returns the total accumulated time for numeric initialization.
Definition at line 256 of file amavanka.hpp.
References FEAT::StopWatch::elapsed().
|
inline |
Returns the total accumulated time for symbolic initialization.
Definition at line 248 of file amavanka.hpp.
References FEAT::StopWatch::elapsed().
|
protected |
deduct macro dofs automatically?
Definition at line 88 of file amavanka.hpp.
Referenced by FEAT::Solver::AmaVanka< Matrix_, Filter_ >::clear_macro_dofs(), and FEAT::Solver::AmaVanka< Matrix_, Filter_ >::push_macro_dofs().
|
protected |
Definition at line 92 of file amavanka.hpp.
|
protected |
the system filter
Definition at line 84 of file amavanka.hpp.
|
protected |
the DOF-macro graphs
Definition at line 92 of file amavanka.hpp.
Referenced by FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::_alloc_device(), FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::_alloc_max_degrees(), FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::_init_numeric_generic(), FEAT::Solver::AmaVanka< Matrix_, Filter_ >::bytes(), FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::init_numeric(), FEAT::Solver::AmaVanka< Matrix_, Filter_ >::init_symbolic(), and FEAT::Solver::AmaVanka< Matrix_, Filter_ >::push_macro_dofs().
|
protected |
the macro mask
Definition at line 94 of file amavanka.hpp.
Referenced by FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::_alloc_device(), and FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::_init_numeric_generic().
|
protected |
the system matrix
Definition at line 82 of file amavanka.hpp.
Referenced by FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::init_numeric().
|
protected |
number of steps
Definition at line 96 of file amavanka.hpp.
Referenced by FEAT::Solver::AmaVanka< Matrix_, Filter_ >::apply().
|
protected |
number of threads for numeric factorization
Definition at line 110 of file amavanka.hpp.
|
protected |
damping parameter
Definition at line 98 of file amavanka.hpp.
Referenced by FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::_init_numeric_cuda(), and FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::_init_numeric_generic().
|
protected |
skip singular macros?
Definition at line 90 of file amavanka.hpp.
Referenced by FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::_alloc_device(), FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::_init_numeric_cuda(), and FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::_init_numeric_generic().
|
protected |
the Vanka preconditioner matrix
Definition at line 86 of file amavanka.hpp.
Referenced by FEAT::Solver::AmaVanka< Matrix_, Filter_ >::bytes(), FEAT::Solver::AmaVanka< Matrix_, Filter_ >::data_size(), and FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::init_numeric().
|
protected |
temporary vectors
Definition at line 100 of file amavanka.hpp.
Referenced by FEAT::Solver::AmaVanka< Matrix_, Filter_ >::bytes().
|
protected |
Definition at line 100 of file amavanka.hpp.
|
protected |
Definition at line 107 of file amavanka.hpp.
|
protected |
Definition at line 105 of file amavanka.hpp.
|
protected |
Definition at line 103 of file amavanka.hpp.