DUNE PDELab (git)

Dune::PDELab::ISTLBackend_OVLP_BCGS_ILUn< GFS, CC > Class Template Reference

Overlapping parallel BiCGStab solver with ILU0 preconditioner. More...

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

Public Member Functions

 ISTLBackend_OVLP_BCGS_ILUn (const GFS &gfs, const CC &cc, int n=1, unsigned maxiter=5000, int verbose=1)
 make a linear solver object More...
 
template<class M , class V , class W >
void apply (M &A, V &z, W &r, typename Dune::template FieldTraits< typename V::ElementType >::real_type reduction)
 solve the given linear system More...
 
template<typename X >
X::ElementType dot (const X &x, const X &y) const
 Dot product of two vectors. It is assumed that the vectors are consistent on the interior+border partition.
 
template<typename X >
Dune::template FieldTraits< typenameX::ElementType >::real_type norm (const X &x) const
 Norm of a right-hand side vector. The vector must be consistent on the interior+border partition.
 
const Dune::PDELab::LinearSolverResult< double > & result () const
 Return access to result data.
 

Detailed Description

template<class GFS, class CC>
class Dune::PDELab::ISTLBackend_OVLP_BCGS_ILUn< GFS, CC >

Overlapping parallel BiCGStab solver with ILU0 preconditioner.

Template Parameters
GFSThe Type of the GridFunctionSpace.
CCThe Type of the Constraints Container.

Constructor & Destructor Documentation

◆ ISTLBackend_OVLP_BCGS_ILUn()

template<class GFS , class CC >
Dune::PDELab::ISTLBackend_OVLP_BCGS_ILUn< GFS, CC >::ISTLBackend_OVLP_BCGS_ILUn ( const GFS &  gfs,
const CC &  cc,
int  n = 1,
unsigned  maxiter = 5000,
int  verbose = 1 
)
inline

make a linear solver object

Parameters
[in]gfsa grid function space
[in]cca constraints container object
[in]nlevel for ILUn
[in]maxitermaximum number of iterations to do
[in]verboseprint messages if true

Member Function Documentation

◆ apply()

template<class GFS , class C , template< class > class Solver>
template<class M , class V , class W >
void Dune::PDELab::ISTLBackend_OVLP_ILUn_Base< GFS, C, Solver >::apply ( M &  A,
V &  z,
W &  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

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 (Nov 23, 23:29, 2024)