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

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

Inheritance diagram for SolverQMRS< VectorType >:
[legend]

Classes

struct  AdditionalData
 
struct  IterationResult
 

Public Member Functions

 SolverQMRS (SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
 
 SolverQMRS (SolverControl &cn, const AdditionalData &data=AdditionalData())
 
template<typename MatrixType , typename PreconditionerType >
void solve (const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
 
virtual void print_vectors (const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
 
- 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

virtual double criterion ()
 

Protected Attributes

VectorType * Vv
 
VectorType * Vx
 
const VectorType * Vb
 
double res2
 
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 , typename PreconditionerType >
IterationResult iterate (const MatrixType &A, const PreconditionerType &precondition)
 

Private Attributes

unsigned int step
 

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 SolverQMRS< VectorType >

Quasi-minimal residual method for symmetric matrices.

The QMRS method is supposed to solve symmetric indefinite linear systems with symmetric, not necessarily definite preconditioners. This version of QMRS is adapted from Freund/Nachtigal: Software for simplified Lanczos and QMR algorithms, Appl. Num. Math. 19 (1995), pp. 319-341

This version is for right preconditioning only, since then only the preconditioner is used: left preconditioning seems to require the inverse.

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.

However, since the QMRS method does not need additional data, the respective structure is empty and does not offer any functionality. The constructor has a default argument, so you may call it without the additional parameter.

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.

Author
Guido Kanschat, 1999

Definition at line 71 of file solver_qmrs.h.

Constructor & Destructor Documentation

◆ SolverQMRS() [1/2]

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

Constructor.

◆ SolverQMRS() [2/2]

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

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

Member Function Documentation

◆ solve()

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

Solve the linear system \(Ax=b\) for x.

◆ print_vectors()

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

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.

◆ criterion()

template<typename VectorType = Vector<double>>
virtual double SolverQMRS< VectorType >::criterion ( )
protectedvirtual

Implementation of the computation of the norm of the residual.

◆ iterate()

template<typename VectorType = Vector<double>>
template<typename MatrixType , typename PreconditionerType >
IterationResult SolverQMRS< 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

◆ Vv

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

Temporary vectors, allocated through the VectorMemory object at the start of the actual solution process and deallocated at the end.

Definition at line 154 of file solver_qmrs.h.

◆ Vx

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

Iteration vector.

Definition at line 162 of file solver_qmrs.h.

◆ Vb

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

RHS vector.

Definition at line 166 of file solver_qmrs.h.

◆ res2

template<typename VectorType = Vector<double>>
double SolverQMRS< VectorType >::res2
protected

Within the iteration loop, the square of the residual vector is stored in this variable. The function criterion uses this variable to compute the convergence value, which in this class is the norm of the residual vector and thus the square root of the res2 value.

Definition at line 174 of file solver_qmrs.h.

◆ additional_data

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

Additional parameters.

Definition at line 179 of file solver_qmrs.h.

◆ step

template<typename VectorType = Vector<double>>
unsigned int SolverQMRS< VectorType >::step
private

Number of the current iteration (accumulated over restarts)

Definition at line 209 of file solver_qmrs.h.


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