DUNE PDELab (git)

Dune::PDELab::ISTLBackend_NOVLP_CG_Jacobi< GFS > Class Template Reference

Nonoverlapping parallel CG solver with Jacobi preconditioner. More...

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

Public Member Functions

 ISTLBackend_NOVLP_CG_Jacobi (const GFS &gfs_, unsigned maxiter_=5000, int verbose_=1)
 make a linear solver object More...
 
template<class V >
V::ElementType norm (const V &v) const
 compute global norm of a vector 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...
 
const LinearSolverResult< double > & result () const
 Return access to result data.
 

Detailed Description

template<class GFS>
class Dune::PDELab::ISTLBackend_NOVLP_CG_Jacobi< GFS >

Nonoverlapping parallel CG solver with Jacobi preconditioner.

Constructor & Destructor Documentation

◆ ISTLBackend_NOVLP_CG_Jacobi()

template<class GFS >
Dune::PDELab::ISTLBackend_NOVLP_CG_Jacobi< GFS >::ISTLBackend_NOVLP_CG_Jacobi ( const GFS &  gfs_,
unsigned  maxiter_ = 5000,
int  verbose_ = 1 
)
inlineexplicit

make a linear solver object

Parameters
gfs_A grid function space
maxiter_Maximum number of iterations to do.
verbose_Verbosity level, directly handed to the CGSolver.

Member Function Documentation

◆ apply()

template<class GFS >
template<class M , class V , class W >
void Dune::PDELab::ISTLBackend_NOVLP_CG_Jacobi< GFS >::apply ( M &  A,
V &  z,
W &  r,
typename Dune::template FieldTraits< typename V::ElementType >::real_type  reduction 
)
inline

solve the given linear system

Parameters
AThe matrix to solve. Should be a matrix from one of PDELabs ISTL backends (only ISTLBCRSMatrixBackend at the moment).
zThe solution vector to be computed
rRight hand side
reductionto be achieved

Solve the linear system A*z=r such that norm(A*z0-r)/norm(A*z-r) < reduction where z0 is the initial value of z.

◆ norm()

template<class GFS >
template<class V >
V::ElementType Dune::PDELab::ISTLBackend_NOVLP_CG_Jacobi< GFS >::norm ( const V &  v) const
inline

compute global norm of a vector

Parameters
vThe vector to compute the norm of. Should be an inconsistent vector (i.e. the entries corresponding a DoF on the border should only contain the summand of this process).

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 12, 23:30, 2024)