DUNE PDELab (2.8)
Dune::PDELab::ISTLBackend_NOVLP_BCGS_NOPREC< GFS > Class Template Reference
Nonoverlapping parallel BiCGStab solver without preconditioner. More...
#include <dune/pdelab/backend/istl/novlpistlsolverbackend.hh>
Public Member Functions | |
ISTLBackend_NOVLP_BCGS_NOPREC (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 Dune::PDELab::LinearSolverResult< double > & | result () const |
Return access to result data. | |
Detailed Description
template<class GFS>
class Dune::PDELab::ISTLBackend_NOVLP_BCGS_NOPREC< GFS >
class Dune::PDELab::ISTLBackend_NOVLP_BCGS_NOPREC< GFS >
Nonoverlapping parallel BiCGStab solver without preconditioner.
Constructor & Destructor Documentation
◆ ISTLBackend_NOVLP_BCGS_NOPREC()
template<class GFS >
|
inlineexplicit |
make a linear solver object
- Parameters
-
[in] gfs_ a grid function space [in] maxiter_ maximum number of iterations to do [in] verbose_ print messages if true
Member Function Documentation
◆ apply()
template<class GFS >
template<class M , class V , class W >
|
inline |
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 GFS >
template<class V >
|
inline |
compute global norm of a vector
- Parameters
-
[in] v the given vector
The documentation for this class was generated from the following file:
- dune/pdelab/backend/istl/novlpistlsolverbackend.hh
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Dec 22, 23:30, 2024)