DUNE PDELab (2.8)

Dune::PDELab::ISTLBackend_SEQ_BCGS_AMG_SSOR< GO > Class Template Reference

Sequential BiCGStab solver preconditioned with AMG smoothed by SSOR. More...

#include <dune/pdelab/backend/istl/seqistlsolverbackend.hh>

Public Member Functions

 ISTLBackend_SEQ_BCGS_AMG_SSOR (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 ISTLAMGStatisticsstatistics () 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_SSOR< GO >

Sequential BiCGStab solver preconditioned with AMG smoothed by SSOR.

Template Parameters
GOThe type of the grid operator (or the fakeGOTraits class for the old grid operator space).

Constructor & Destructor Documentation

◆ ISTLBackend_SEQ_BCGS_AMG_SSOR()

template<class GO >
Dune::PDELab::ISTLBackend_SEQ_BCGS_AMG_SSOR< GO >::ISTLBackend_SEQ_BCGS_AMG_SSOR ( unsigned  maxiter_ = 5000,
int  verbose_ = 1,
bool  reuse_ = false,
bool  usesuperlu_ = true 
)
inline

Constructor.

Parameters
maxiter_The maximum number of iterations allowed.
verbose_The verbosity level to use.
reuse_Set true, if the Matrix to be used is always identical (AMG aggregation is then only performed once).
usesuperlu_Set false, to suppress the no SuperLU warning

Member Function Documentation

◆ apply()

template<class GO , template< class, class, class, int > class Preconditioner, template< class > class Solver, bool skipBlocksizeCheck = false>
void Dune::PDELab::ISTLBackend_SEQ_AMG< GO, Preconditioner, Solver, skipBlocksizeCheck >::apply ( M &  A,
V &  z,
V &  r,
typename Dune::template FieldTraits< typename V::ElementType >::real_type  reduction 
)
inlineinherited

solve the given linear system

Parameters
[in]Athe given matrix
[out]zthe solution vector to be computed
[in]rright hand side
[in]reductionto be achieved

◆ norm()

template<class GO , template< class, class, class, int > class Preconditioner, template< class > class Solver, bool skipBlocksizeCheck = false>
V::ElementType Dune::PDELab::ISTLBackend_SEQ_AMG< GO, Preconditioner, Solver, skipBlocksizeCheck >::norm ( const V &  v) const
inlineinherited

compute global norm of a vector

Parameters
[in]vthe given vector

◆ setparams()

template<class GO , template< class, class, class, int > class Preconditioner, template< class > class Solver, bool skipBlocksizeCheck = false>
void Dune::PDELab::ISTLBackend_SEQ_AMG< GO, Preconditioner, Solver, skipBlocksizeCheck >::setparams ( Parameters  params_)
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>
const ISTLAMGStatistics & Dune::PDELab::ISTLBackend_SEQ_AMG< GO, Preconditioner, Solver, skipBlocksizeCheck >::statistics ( ) const
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:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)