DUNE PDELab (git)

Dune::PDELab::ISTLBackend_SEQ_BCGS_ILU0 Class Reference

Backend for sequential BiCGSTAB solver with ILU0 preconditioner. More...

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

Public Member Functions

 ISTLBackend_SEQ_BCGS_ILU0 (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 W::ElementType >::real_type reduction)
 solve the given linear system More...
 
template<class V >
Dune::template FieldTraits< typenameV::ElementType >::real_type norm (const V &v) const
 compute global norm of a vector More...
 
const Dune::PDELab::LinearSolverResult< double > & result () const
 Return access to result data.
 

Detailed Description

Backend for sequential BiCGSTAB solver with ILU0 preconditioner.

Constructor & Destructor Documentation

◆ ISTLBackend_SEQ_BCGS_ILU0()

Dune::PDELab::ISTLBackend_SEQ_BCGS_ILU0::ISTLBackend_SEQ_BCGS_ILU0 ( unsigned  maxiter_ = 5000,
int  verbose_ = 1 
)
inlineexplicit

make a linear solver object

Parameters
[in]maxiter_maximum number of iterations to do
[in]verbose_print messages if true

Member Function Documentation

◆ apply()

template<template< typename > class Solver>
template<class M , class V , class W >
void Dune::PDELab::ISTLBackend_SEQ_ILU0< Solver >::apply ( M &  A,
V &  z,
W &  r,
typename Dune::template FieldTraits< typename W::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 V >
Dune::template FieldTraits< typenameV::ElementType >::real_type Dune::PDELab::SequentialNorm::norm ( const V &  v) const
inlineinherited

compute global norm of a vector

\param[in] v the given vector

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)