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

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

Inheritance diagram for IterativeInverse< VectorType >:
[legend]

Public Member Functions

template<typename MatrixType , typename PreconditionerType >
void initialize (const MatrixType &, const PreconditionerType &)
 
void clear ()
 
void vmult (VectorType &dst, const VectorType &src) const
 
template<class OtherVectorType >
void vmult (OtherVectorType &dst, const OtherVectorType &src) const
 
- 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)
 

Public Attributes

SolverSelector< VectorType > solver
 

Private Attributes

std_cxx11::shared_ptr< PointerMatrixBase< VectorType > > matrix
 
std_cxx11::shared_ptr< PointerMatrixBase< VectorType > > preconditioner
 

Additional Inherited Members

- 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>
class IterativeInverse< VectorType >

Implementation of the inverse of a matrix, using an iterative method.

The function vmult() of this class starts an iterative solver in order to approximate the action of the inverse matrix.

Krylov space methods like SolverCG or SolverBicgstab become inefficient if solution down to machine accuracy is needed. This is due to the fact, that round-off errors spoil the orthogonality of the vector sequences. Therefore, a nested iteration of two methods is proposed: The outer method is SolverRichardson, since it is robust with respect to round-of errors. The inner loop is an appropriate Krylov space method, since it is fast.

// Declare related objects
ReductionControl inner_control (10, 1.e-30, 1.e-2);
inner_precondition.initialize (A, 1.2);
precondition.initialize (A, inner_precondition);
precondition.solver.select("cg");
precondition.solver.set_control(inner_control);
SolverControl outer_control(100, 1.e-16);
outer_iteration.solve (A, x, b, precondition);

Each time we call the inner loop, reduction of the residual by a factor 1.e-2 is attempted. Since the right hand side vector of the inner iteration is the residual of the outer loop, the relative errors are far from machine accuracy, even if the errors of the outer loop are in the range of machine accuracy.

Deprecated:
If deal.II was configured with C++11 support, use the LinearOperator class instead, see the module on linear operators for further details.
Author
Guido Kanschat
Date
2010

Definition at line 79 of file iterative_inverse.h.

Member Function Documentation

◆ initialize()

template<typename VectorType >
template<typename MatrixType , typename PreconditionerType >
void IterativeInverse< VectorType >::initialize ( const MatrixType &  m,
const PreconditionerType &  p 
)
inline

Initialization function. Provide a matrix and preconditioner for the solve in vmult().

Definition at line 130 of file iterative_inverse.h.

◆ clear()

template<typename VectorType >
void IterativeInverse< VectorType >::clear ( )
inline

Delete the pointers to matrix and preconditioner.

Definition at line 141 of file iterative_inverse.h.

◆ vmult() [1/2]

template<typename VectorType >
void IterativeInverse< VectorType >::vmult ( VectorType &  dst,
const VectorType &  src 
) const
inline

Solve for right hand side src.

Definition at line 150 of file iterative_inverse.h.

◆ vmult() [2/2]

template<typename VectorType >
template<class OtherVectorType >
void IterativeInverse< VectorType >::vmult ( OtherVectorType &  dst,
const OtherVectorType &  src 
) const
inline

Solve for right hand side src, but allow for the fact that the vectors given to this function have different type from the vectors used by the inner solver.

Definition at line 162 of file iterative_inverse.h.

Member Data Documentation

◆ solver

template<typename VectorType>
SolverSelector<VectorType> IterativeInverse< VectorType >::solver

The solver, which allows selection of the actual solver as well as adjustment of parameters.

Definition at line 111 of file iterative_inverse.h.

◆ matrix

template<typename VectorType>
std_cxx11::shared_ptr<PointerMatrixBase<VectorType> > IterativeInverse< VectorType >::matrix
private

The matrix in use.

Definition at line 117 of file iterative_inverse.h.

◆ preconditioner

template<typename VectorType>
std_cxx11::shared_ptr<PointerMatrixBase<VectorType> > IterativeInverse< VectorType >::preconditioner
private

The preconditioner to use.

Definition at line 122 of file iterative_inverse.h.


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