Reference documentation for deal.II version 8.5.1
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
SolverBicgstab< VectorType > Class Template Reference

#include <deal.II/lac/solver_bicgstab.h>

Inheritance diagram for SolverBicgstab< VectorType >:
[legend]

Classes

struct  AdditionalData
 
struct  IterationResult
 

Public Member Functions

 SolverBicgstab (SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
 
 SolverBicgstab (SolverControl &cn, const AdditionalData &data=AdditionalData())
 
virtual ~SolverBicgstab ()
 
template<typename MatrixType , typename PreconditionerType >
void solve (const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
 
- Public Member Functions inherited from Solver< VectorType >
 Solver (SolverControl &solver_control, VectorMemory< VectorType > &vector_memory)
 
 Solver (SolverControl &solver_control)
 
boost::signals2::connection connect (const std_cxx11::function< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType &current_iterate)> &slot)
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&)
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&)
 
void subscribe (const char *identifier=0) const
 
void unsubscribe (const char *identifier=0) const
 
unsigned int n_subscriptions () const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Protected Member Functions

template<typename MatrixType >
double criterion (const MatrixType &A, const VectorType &x, const VectorType &b)
 
virtual void print_vectors (const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
 

Protected Attributes

VectorType * Vx
 
VectorType * Vr
 
VectorType * Vrbar
 
VectorType * Vp
 
VectorType * Vy
 
VectorType * Vz
 
VectorType * Vt
 
VectorType * Vv
 
const VectorType * Vb
 
double alpha
 
double beta
 
double omega
 
double rho
 
double rhobar
 
unsigned int step
 
double res
 
AdditionalData additional_data
 
- Protected Attributes inherited from Solver< VectorType >
GrowingVectorMemory< VectorType > static_vector_memory
 
VectorMemory< VectorType > & memory
 
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType &current_iterate), StateCombineriteration_status
 

Private Member Functions

template<typename MatrixType >
SolverControl::State start (const MatrixType &A)
 
template<typename MatrixType , typename PreconditionerType >
IterationResult iterate (const MatrixType &A, const PreconditionerType &precondition)
 

Additional Inherited Members

- Public Types inherited from Solver< VectorType >
typedef VectorType vector_type
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, char *arg2, std::string &arg3)
 
static ::ExceptionBaseExcNoSubscriber (char *arg1, char *arg2)
 

Detailed Description

template<typename VectorType = Vector<double>>
class SolverBicgstab< VectorType >

Bicgstab algorithm by van der Vorst.

For the requirements on matrices and vectors in order to work with this class, see the documentation of the Solver base class.

Like all other solver classes, this class has a local structure called AdditionalData which is used to pass additional parameters to the solver, like damping parameters or the number of temporary vectors. We use this additional structure instead of passing these values directly to the constructor because this makes the use of the SolverSelector and other classes much easier and guarantees that these will continue to work even if number or type of the additional parameters for a certain solver changes.

The Bicgstab-method has two additional parameters: the first is a boolean, deciding whether to compute the actual residual in each step (true) or to use the length of the computed orthogonal residual (false). Note that computing the residual causes a third matrix-vector-multiplication, though no additional preconditioning, in each step. The reason for doing this is, that the size of the orthogonalized residual computed during the iteration may be larger by orders of magnitude than the true residual. This is due to numerical instabilities related to badly conditioned matrices. Since this instability results in a bad stopping criterion, the default for this parameter is true. Whenever the user knows that the estimated residual works reasonably as well, the flag should be set to false in order to increase the performance of the solver.

The second parameter is the size of a breakdown criterion. It is difficult to find a general good criterion, so if things do not work for you, try to change this value.

Observing the progress of linear solver iterations

The solve() function of this class uses the mechanism described in the Solver base class to determine convergence. This mechanism can also be used to observe the progress of the iteration.

Definition at line 73 of file solver_bicgstab.h.

Constructor & Destructor Documentation

◆ SolverBicgstab() [1/2]

template<typename VectorType = Vector<double>>
SolverBicgstab< VectorType >::SolverBicgstab ( SolverControl cn,
VectorMemory< VectorType > &  mem,
const AdditionalData data = AdditionalData() 
)

Constructor.

◆ SolverBicgstab() [2/2]

template<typename VectorType = Vector<double>>
SolverBicgstab< VectorType >::SolverBicgstab ( SolverControl cn,
const AdditionalData data = AdditionalData() 
)

Constructor. Use an object of type GrowingVectorMemory as a default to allocate memory.

◆ ~SolverBicgstab()

template<typename VectorType = Vector<double>>
virtual SolverBicgstab< VectorType >::~SolverBicgstab ( )
virtual

Virtual destructor.

Member Function Documentation

◆ solve()

template<typename VectorType = Vector<double>>
template<typename MatrixType , typename PreconditionerType >
void SolverBicgstab< VectorType >::solve ( const MatrixType &  A,
VectorType &  x,
const VectorType &  b,
const PreconditionerType &  precondition 
)

Solve primal problem only.

◆ criterion()

template<typename VectorType = Vector<double>>
template<typename MatrixType >
double SolverBicgstab< VectorType >::criterion ( const MatrixType &  A,
const VectorType &  x,
const VectorType &  b 
)
protected

Computation of the stopping criterion.

◆ print_vectors()

template<typename VectorType = Vector<double>>
virtual void SolverBicgstab< VectorType >::print_vectors ( const unsigned int  step,
const VectorType &  x,
const VectorType &  r,
const VectorType &  d 
) const
protectedvirtual

Interface for derived class. This function gets the current iteration vector, the residual and the update vector in each step. It can be used for a graphical output of the convergence history.

◆ start()

template<typename VectorType = Vector<double>>
template<typename MatrixType >
SolverControl::State SolverBicgstab< VectorType >::start ( const MatrixType &  A)
private

Everything before the iteration loop.

◆ iterate()

template<typename VectorType = Vector<double>>
template<typename MatrixType , typename PreconditionerType >
IterationResult SolverBicgstab< VectorType >::iterate ( const MatrixType &  A,
const PreconditionerType &  precondition 
)
private

The iteration loop itself. The function returns a structure indicating what happened in this function.

Member Data Documentation

◆ Vx

template<typename VectorType = Vector<double>>
VectorType* SolverBicgstab< VectorType >::Vx
protected

Auxiliary vector.

Definition at line 159 of file solver_bicgstab.h.

◆ Vr

template<typename VectorType = Vector<double>>
VectorType* SolverBicgstab< VectorType >::Vr
protected

Auxiliary vector.

Definition at line 163 of file solver_bicgstab.h.

◆ Vrbar

template<typename VectorType = Vector<double>>
VectorType* SolverBicgstab< VectorType >::Vrbar
protected

Auxiliary vector.

Definition at line 167 of file solver_bicgstab.h.

◆ Vp

template<typename VectorType = Vector<double>>
VectorType* SolverBicgstab< VectorType >::Vp
protected

Auxiliary vector.

Definition at line 171 of file solver_bicgstab.h.

◆ Vy

template<typename VectorType = Vector<double>>
VectorType* SolverBicgstab< VectorType >::Vy
protected

Auxiliary vector.

Definition at line 175 of file solver_bicgstab.h.

◆ Vz

template<typename VectorType = Vector<double>>
VectorType* SolverBicgstab< VectorType >::Vz
protected

Auxiliary vector.

Definition at line 179 of file solver_bicgstab.h.

◆ Vt

template<typename VectorType = Vector<double>>
VectorType* SolverBicgstab< VectorType >::Vt
protected

Auxiliary vector.

Definition at line 183 of file solver_bicgstab.h.

◆ Vv

template<typename VectorType = Vector<double>>
VectorType* SolverBicgstab< VectorType >::Vv
protected

Auxiliary vector.

Definition at line 187 of file solver_bicgstab.h.

◆ Vb

template<typename VectorType = Vector<double>>
const VectorType* SolverBicgstab< VectorType >::Vb
protected

Right hand side vector.

Definition at line 191 of file solver_bicgstab.h.

◆ alpha

template<typename VectorType = Vector<double>>
double SolverBicgstab< VectorType >::alpha
protected

Auxiliary value.

Definition at line 196 of file solver_bicgstab.h.

◆ beta

template<typename VectorType = Vector<double>>
double SolverBicgstab< VectorType >::beta
protected

Auxiliary value.

Definition at line 200 of file solver_bicgstab.h.

◆ omega

template<typename VectorType = Vector<double>>
double SolverBicgstab< VectorType >::omega
protected

Auxiliary value.

Definition at line 204 of file solver_bicgstab.h.

◆ rho

template<typename VectorType = Vector<double>>
double SolverBicgstab< VectorType >::rho
protected

Auxiliary value.

Definition at line 208 of file solver_bicgstab.h.

◆ rhobar

template<typename VectorType = Vector<double>>
double SolverBicgstab< VectorType >::rhobar
protected

Auxiliary value.

Definition at line 212 of file solver_bicgstab.h.

◆ step

template<typename VectorType = Vector<double>>
unsigned int SolverBicgstab< VectorType >::step
protected

Current iteration step.

Definition at line 217 of file solver_bicgstab.h.

◆ res

template<typename VectorType = Vector<double>>
double SolverBicgstab< VectorType >::res
protected

Residual.

Definition at line 222 of file solver_bicgstab.h.

◆ additional_data

template<typename VectorType = Vector<double>>
AdditionalData SolverBicgstab< VectorType >::additional_data
protected

Additional parameters.

Definition at line 227 of file solver_bicgstab.h.


The documentation for this class was generated from the following file: