DUNE PDELab (git)
Dune::PDELab::ISTLBackend_SEQ_BCGS_AMG_SOR< GO > Class Template Reference
Sequential BiCGSTAB solver preconditioned with AMG smoothed by SOR. More...
#include <dune/pdelab/backend/istl/seqistlsolverbackend.hh>
Public Member Functions | |
ISTLBackend_SEQ_BCGS_AMG_SOR (unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true) | |
Constructor. More... | |
void | setparams (Parameters params_) |
set AMG parameters More... | |
void | setReuse (bool reuse_) |
Set whether the AMG should be reused again during call to apply(). | |
bool | getReuse () const |
Return whether the AMG is reused during call to apply() | |
V::ElementType | norm (const V &v) const |
compute global norm of a vector More... | |
void | apply (M &A, V &z, V &r, typename Dune::template FieldTraits< typename V::ElementType >::real_type reduction) |
solve the given linear system More... | |
const ISTLAMGStatistics & | statistics () const |
Get statistics of the AMG solver (no of levels, timings). More... | |
const Dune::PDELab::LinearSolverResult< double > & | result () const |
Return access to result data. | |
Detailed Description
template<class GO>
class Dune::PDELab::ISTLBackend_SEQ_BCGS_AMG_SOR< GO >
class Dune::PDELab::ISTLBackend_SEQ_BCGS_AMG_SOR< GO >
Sequential BiCGSTAB solver preconditioned with AMG smoothed by SOR.
- Template Parameters
-
GO The type of the grid operator (or the fakeGOTraits class for the old grid operator space).
Constructor & Destructor Documentation
◆ ISTLBackend_SEQ_BCGS_AMG_SOR()
template<class GO >
|
inline |
Member Function Documentation
◆ apply()
template<class GO , template< class, class, class, int > class Preconditioner, template< class > class Solver, bool skipBlocksizeCheck = false>
|
inlineinherited |
solve the given linear system
- Parameters
-
[in] A the given matrix [out] z the solution vector to be computed [in] r right hand side [in] reduction to be achieved
◆ norm()
template<class GO , template< class, class, class, int > class Preconditioner, template< class > class Solver, bool skipBlocksizeCheck = false>
|
inlineinherited |
compute global norm of a vector
- Parameters
-
[in] v the given vector
◆ setparams()
template<class GO , template< class, class, class, int > class Preconditioner, template< class > class Solver, bool skipBlocksizeCheck = false>
|
inlineinherited |
set AMG parameters
- Parameters
-
[in] params_ a parameter object of Type Dune::Amg::Parameters
◆ statistics()
template<class GO , template< class, class, class, int > class Preconditioner, template< class > class Solver, bool skipBlocksizeCheck = false>
|
inlineinherited |
Get statistics of the AMG solver (no of levels, timings).
- Returns
- statistis of the AMG solver.
The documentation for this class was generated from the following file:
- dune/pdelab/backend/istl/seqistlsolverbackend.hh
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Nov 12, 23:30, 2024)