|
FEAT 3
Finite Element Analysis Toolbox
|
Multigrid preconditioner implementation. More...
#include <multigrid.hpp>
Public Types | |
| typedef SolverBase< typename SystemMatrix_::VectorTypeR > | BaseClass |
| our base-class More... | |
| typedef MatrixType::DataType | DataType |
| our data type More... | |
| typedef LevelType::SystemFilterType | FilterType |
| the system filter type More... | |
| typedef MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ > | HierarchyType |
| our compatible multigrid hierarchy class More... | |
| typedef HierarchyType::LevelInfo | LevelInfo |
| the level info type More... | |
| typedef HierarchyType::LevelType | LevelType |
| the level type More... | |
| typedef LevelType::SystemMatrixType | MatrixType |
| the system matrix type More... | |
| typedef LevelType::SolverType | SolverType |
| the sub-solver type More... | |
| typedef LevelType::TransferOperatorType | TransferOperatorType |
| the transfer operator type More... | |
| typedef LevelType::SystemVectorType | VectorType |
| the system vector type More... | |
Public Member Functions | |
| MultiGrid (std::shared_ptr< HierarchyType > hierarchy, MultiGridCycle cycle, int top_level=0, int crs_level=-1) | |
| Constructor. More... | |
| virtual | ~MultiGrid () |
| virtual destructor More... | |
| virtual Status | apply (VectorType &vec_cor, const VectorType &vec_def) override |
| Applies the multigrid preconditioner. More... | |
| virtual void | done () |
| Finalization method. More... | |
| virtual void | done_numeric () override |
| Numeric finalization function. More... | |
| virtual void | done_symbolic () override |
| Symbolic finalization function. More... | |
| Index | get_crs_level () const |
| Returns the currently selected coarse-level index. More... | |
| MultiGridCycle | get_cycle () const |
| Returns the currently selected cycle. More... | |
| Index | get_top_level () const |
| Returns the currently selected top-level index. More... | |
| virtual void | init () |
| Initialization method. More... | |
| virtual void | init_numeric () override |
| Numeric initialization function. More... | |
| virtual void | init_symbolic () override |
| Symbolic initialization function. More... | |
| virtual String | name () const override |
| Returns a descriptive string. More... | |
| void | set_adapt_cgc (MultiGridAdaptCGC adapt_cgc) |
| Sets the adaption mode for the coarse grid correction. More... | |
| void | set_cycle (MultiGridCycle cycle) |
| Sets a new cycle. More... | |
| void | set_levels (int top_level, int crs_level) |
| Sets the level range for this multigrid. More... | |
Protected Member Functions | |
| Status | _apply_coarse () |
| Applies the coarse grid solver on the coarse level. More... | |
| Status | _apply_cycle_f () |
| Applies a single F-cycle. More... | |
| Status | _apply_cycle_v () |
| Applies a single V-cycle. More... | |
| Status | _apply_cycle_w () |
| Applies a single W-cycle. More... | |
| virtual Status | _apply_prol (const Index cur_lvl, bool cur_smooth) |
| Prolongates from the coarse level onto the current level. More... | |
| Status | _apply_rest (const Index cur_lvl, bool cur_smooth) |
| Restricts from the current level onto the coarse level. More... | |
| bool | _apply_smooth_def (Index cur_lvl, SolverType &smoother) |
| Status | _apply_smooth_peak (Index cur_lvl) |
| Applies the peak-smoother on the current level. More... | |
Protected Attributes | |
| MultiGridAdaptCGC | _adapt_cgc |
| the coarse grid correction type More... | |
| std::vector< int > | _counters |
| W-cycle level counter vector. More... | |
| Index | _crs_level |
| the coarse-level of this multigrid More... | |
| MultiGridCycle | _cycle |
| the multigrid cycle More... | |
| std::shared_ptr< HierarchyType > | _hierarchy |
| the multigrid hierarchy object More... | |
| std::vector< double > | _mpi_execs_blas2 |
| array containing toe of mpi execution blas2 for each processed level More... | |
| std::vector< double > | _mpi_execs_blas3 |
| array containing toe of mpi execution blas3 for each processed level More... | |
| std::vector< double > | _mpi_execs_collective |
| array containing toe of mpi execution collective for each processed level More... | |
| std::vector< double > | _mpi_execs_reduction |
| array containing toe of mpi execution reduction for each processed level More... | |
| std::vector< double > | _mpi_waits_blas2 |
| array containing toe of mpi wait blas2 for each processed level More... | |
| std::vector< double > | _mpi_waits_blas3 |
| array containing toe of mpi wait blas3 for each processed level More... | |
| std::vector< double > | _mpi_waits_collective |
| array containing toe of mpi wait collective for each processed level More... | |
| std::vector< double > | _mpi_waits_reduction |
| array containing toe of mpi wait reduction for each processed level More... | |
| std::vector< double > | _toes |
| array containing toe for each processed level More... | |
| Index | _top_level |
| the top-level of this multigrid More... | |
Multigrid preconditioner implementation.
| SystemMatrix_ | The class representing the system matrix. |
| SystemFilter_ | The class representing the system filter. |
| ProlOperator_ | The class representing the prolongation operator (usually a matrix type). |
| RestOperator_ | The class representing the restriction operator (usually a matrix type). |
Definition at line 982 of file multigrid.hpp.
| typedef SolverBase<typename SystemMatrix_::VectorTypeR> FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::BaseClass |
our base-class
Definition at line 987 of file multigrid.hpp.
| typedef MatrixType::DataType FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::DataType |
our data type
Definition at line 1009 of file multigrid.hpp.
| typedef LevelType::SystemFilterType FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::FilterType |
the system filter type
Definition at line 1000 of file multigrid.hpp.
| typedef MultiGridHierarchy<SystemMatrix_, SystemFilter_, TransferOperator_> FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::HierarchyType |
our compatible multigrid hierarchy class
Definition at line 990 of file multigrid.hpp.
| typedef HierarchyType::LevelInfo FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo |
the level info type
Definition at line 995 of file multigrid.hpp.
| typedef HierarchyType::LevelType FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelType |
the level type
Definition at line 993 of file multigrid.hpp.
| typedef LevelType::SystemMatrixType FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::MatrixType |
the system matrix type
Definition at line 998 of file multigrid.hpp.
| typedef LevelType::SolverType FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::SolverType |
the sub-solver type
Definition at line 1006 of file multigrid.hpp.
| typedef LevelType::TransferOperatorType FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::TransferOperatorType |
the transfer operator type
Definition at line 1004 of file multigrid.hpp.
| typedef LevelType::SystemVectorType FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::VectorType |
the system vector type
Definition at line 1002 of file multigrid.hpp.
|
inlineexplicit |
Constructor.
| [in] | hierarchy | A pointer to the multigrid hierarchy object. |
| [in] | cycle | The desired multigrid cycle. |
| [in] | top_level | The desired top-level for this multigrid solver. Set to 0 to use the finest level in the multigrid hierarchy. |
| [in] | crs_level | The desired coarse-level for this multigrid solver. Set to -1 to use the coarsest level in the multigrid hierarchy. |
Definition at line 1066 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_crs_level, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_top_level, and XASSERTM.
|
inlinevirtual |
virtual destructor
Definition at line 1079 of file multigrid.hpp.
|
inlineprotected |
Applies the coarse grid solver on the coarse level.
Definition at line 1524 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_crs_level, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_hierarchy, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_toes, FEAT::Solver::aborted, FEAT::TimeStamp::elapsed_now(), FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::level, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::name(), FEAT::Solver::status_success(), FEAT::Solver::success, FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::time_coarse, and FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::vec_rhs.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_f(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_v(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_w().
|
inlineprotected |
Applies a single F-cycle.
Definition at line 1354 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_smooth_peak(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_crs_level, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_hierarchy, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_top_level, FEAT::Math::min(), and FEAT::Solver::success.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply().
|
inlineprotected |
Applies a single V-cycle.
Definition at line 1333 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), and FEAT::Solver::success.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply().
|
inlineprotected |
Applies a single W-cycle.
Definition at line 1410 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_smooth_peak(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_counters, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_crs_level, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_hierarchy, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_top_level, FEAT::Math::min(), FEAT::Solver::success, and XASSERTM.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply().
|
inlineprotectedvirtual |
Prolongates from the coarse level onto the current level.
| [in] | cur_lvl | The level onto which to prolongate. |
| [in] | cur_smooth | Specifies whether to apply the post-smoother on the current level. |
Definition at line 1843 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_adapt_cgc, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_crs_level, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_hierarchy, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_toes, FEAT::TimeStamp::elapsed_now(), FEAT::Solver::Fixed, FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::level, FEAT::Math::min(), FEAT::Solver::MinDefect, FEAT::Solver::MinEnergy, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::name(), FEAT::Solver::success, FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::time_defect, FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::time_smooth, FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::time_transfer, FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::vec_rhs, and XASSERTM.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_f(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_v(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_w().
|
inlineprotected |
Restricts from the current level onto the coarse level.
| [in] | cur_lvl | The level from which to restrict. |
| [in] | cur_smooth | Specifies whether to apply the pre-smoother on the current level. |
Definition at line 1711 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_crs_level, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_hierarchy, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_toes, FEAT::TimeStamp::elapsed_now(), FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::level, FEAT::Math::min(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::name(), FEAT::Solver::success, FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::time_defect, FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::time_smooth, FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::time_transfer, FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::vec_rhs, and XASSERTM.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_f(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_v(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_w().
|
inlineprotected |
Definition at line 1589 of file multigrid.hpp.
|
inlineprotected |
Applies the peak-smoother on the current level.
| [in] | cur_lvl | The current level. |
Definition at line 1630 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_hierarchy, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_toes, FEAT::Solver::aborted, FEAT::TimeStamp::elapsed(), FEAT::TimeStamp::elapsed_now(), FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::level, FEAT::Solver::success, FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::time_defect, and FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::vec_rhs.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_f(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_w().
|
inlineoverridevirtual |
Applies the multigrid preconditioner.
| [out] | vec_cor | The vector that shall receive the solution of the linear system. It is assumed to be allocated, but its numerical contents may be undefined upon calling this method. |
| [in] | vec_def | The vector that represents the right-hand-side of the linear system to be solved. |
Implements FEAT::Solver::SolverBase< SystemMatrix_::VectorTypeR >.
Definition at line 1268 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_f(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_v(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_w(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_cycle, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_hierarchy, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_toes, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_top_level, FEAT::Solver::aborted, FEAT::Solver::F, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::name(), FEAT::Solver::success, FEAT::Solver::undefined, FEAT::Solver::V, FEAT::Solver::MultiGridHierarchy< SystemMatrix_, SystemFilter_, TransferOperator_ >::LevelInfo::vec_rhs, and FEAT::Solver::W.
|
inlinevirtualinherited |
|
inlineoverridevirtual |
Numeric finalization function.
Reimplemented from FEAT::Solver::SolverBase< SystemMatrix_::VectorTypeR >.
Definition at line 1250 of file multigrid.hpp.
References FEAT::Solver::SolverBase< Vector_ >::done_numeric().
|
inlineoverridevirtual |
Symbolic finalization function.
Reimplemented from FEAT::Solver::SolverBase< SystemMatrix_::VectorTypeR >.
Definition at line 1208 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_counters, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_toes, and FEAT::Solver::SolverBase< Vector_ >::done_symbolic().
|
inline |
Returns the currently selected coarse-level index.
Definition at line 1150 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_crs_level.
|
inline |
Returns the currently selected cycle.
Definition at line 1110 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_cycle.
|
inline |
Returns the currently selected top-level index.
Definition at line 1140 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_top_level.
|
inlinevirtualinherited |
|
inlineoverridevirtual |
Numeric initialization function.
Reimplemented from FEAT::Solver::SolverBase< SystemMatrix_::VectorTypeR >.
Definition at line 1233 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_hierarchy, FEAT::Solver::SolverBase< Vector_ >::init_numeric(), and XASSERTM.
|
inlineoverridevirtual |
Symbolic initialization function.
Reimplemented from FEAT::Solver::SolverBase< SystemMatrix_::VectorTypeR >.
Definition at line 1174 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_counters, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_hierarchy, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_execs_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas2, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_blas3, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_collective, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_mpi_waits_reduction, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_toes, FEAT::Solver::SolverBase< Vector_ >::init_symbolic(), and XASSERTM.
|
inlineoverridevirtual |
Returns a descriptive string.
Implements FEAT::Solver::SolverBase< SystemMatrix_::VectorTypeR >.
Definition at line 1160 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_cycle, and FEAT::stringify().
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply().
|
inline |
Sets the adaption mode for the coarse grid correction.
| [in] | adapt_cgc | The adaption mode for the coarse grid correction. |
Definition at line 1089 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_adapt_cgc.
|
inline |
Sets a new cycle.
| [in] | cycle | The new cycle to be used. |
Definition at line 1100 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_cycle.
|
inline |
Sets the level range for this multigrid.
| [in] | top_level | The desired top-level for this multigrid solver. Set to 0 to use the finest level in the multigrid hierarchy. |
| [in] | crs_level | The desired coarse-level for this multigrid solver. Set to -1 to use the coarsest level in the multigrid hierarchy. |
Definition at line 1126 of file multigrid.hpp.
References FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_crs_level, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_hierarchy, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_top_level, and XASSERTM.
|
protected |
the coarse grid correction type
Definition at line 1017 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::set_adapt_cgc().
|
protected |
W-cycle level counter vector.
Definition at line 1024 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_w(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::done_symbolic(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_symbolic().
|
protected |
the coarse-level of this multigrid
Definition at line 1021 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::MultiGrid(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_f(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_w(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::get_crs_level(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::set_levels().
|
protected |
the multigrid cycle
Definition at line 1015 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::get_cycle(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::name(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::set_cycle().
|
protected |
the multigrid hierarchy object
Definition at line 1013 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_f(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_w(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_smooth_peak(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_numeric(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_symbolic(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::set_levels().
|
protected |
array containing toe of mpi execution blas2 for each processed level
Definition at line 1031 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_smooth_peak(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::done_symbolic(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_symbolic().
|
protected |
array containing toe of mpi execution blas3 for each processed level
Definition at line 1033 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_smooth_peak(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::done_symbolic(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_symbolic().
|
protected |
array containing toe of mpi execution collective for each processed level
Definition at line 1035 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_smooth_peak(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::done_symbolic(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_symbolic().
|
protected |
array containing toe of mpi execution reduction for each processed level
Definition at line 1029 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_smooth_peak(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::done_symbolic(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_symbolic().
|
protected |
array containing toe of mpi wait blas2 for each processed level
Definition at line 1039 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_smooth_peak(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::done_symbolic(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_symbolic().
|
protected |
array containing toe of mpi wait blas3 for each processed level
Definition at line 1041 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_smooth_peak(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::done_symbolic(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_symbolic().
|
protected |
array containing toe of mpi wait collective for each processed level
Definition at line 1043 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_smooth_peak(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::done_symbolic(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_symbolic().
|
protected |
array containing toe of mpi wait reduction for each processed level
Definition at line 1037 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_smooth_peak(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::done_symbolic(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_symbolic().
|
protected |
array containing toe for each processed level
Definition at line 1027 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_coarse(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_prol(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_rest(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_smooth_peak(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::done_symbolic(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_symbolic().
|
protected |
the top-level of this multigrid
Definition at line 1019 of file multigrid.hpp.
Referenced by FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::MultiGrid(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_f(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::_apply_cycle_w(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::apply(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::get_top_level(), and FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::set_levels().