Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Classes | Public Member Functions | Static Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
SolverGMRES< VectorType > Class Template Reference

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

Inheritance diagram for SolverGMRES< VectorType >:
[legend]

Classes

struct  AdditionalData
 

Public Member Functions

 SolverGMRES (SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
 
 SolverGMRES (SolverControl &cn, const AdditionalData &data=AdditionalData())
 
 SolverGMRES (const SolverGMRES< VectorType > &)=delete
 
template<typename MatrixType , typename PreconditionerType >
void solve (const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
 
boost::signals2::connection connect_condition_number_slot (const std::function< void(double)> &slot, const bool every_iteration=false)
 
boost::signals2::connection connect_eigenvalues_slot (const std::function< void(const std::vector< std::complex< double >> &)> &slot, const bool every_iteration=false)
 
boost::signals2::connection connect_hessenberg_slot (const std::function< void(const FullMatrix< double > &)> &slot, const bool every_iteration=true)
 
boost::signals2::connection connect_krylov_space_slot (const std::function< void(const internal::SolverGMRESImplementation::TmpVectors< VectorType > &)> &slot)
 
boost::signals2::connection connect_re_orthogonalization_slot (const std::function< void(int)> &slot)
 
- Public Member Functions inherited from SolverBase< Vector< double > >
 SolverBase (SolverControl &solver_control, VectorMemory< Vector< double > > &vector_memory)
 
 SolverBase (SolverControl &solver_control)
 
boost::signals2::connection connect (const std::function< SolverControl::State(const unsigned int iteration, const double check_value, const Vector< double > &current_iterate)> &slot)
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Static Public Member Functions

static ::ExceptionBaseExcTooFewTmpVectors (int arg1)
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Protected Member Functions

virtual double criterion ()
 
void givens_rotation (Vector< double > &h, Vector< double > &b, Vector< double > &ci, Vector< double > &si, int col) const
 

Static Protected Member Functions

static double modified_gram_schmidt (const internal::SolverGMRESImplementation::TmpVectors< VectorType > &orthogonal_vectors, const unsigned int dim, const unsigned int accumulated_iterations, VectorType &vv, Vector< double > &h, bool &re_orthogonalize, const boost::signals2::signal< void(int)> &re_orthogonalize_signal=boost::signals2::signal< void(int)>())
 
static void compute_eigs_and_cond (const FullMatrix< double > &H_orig, const unsigned int dim, const boost::signals2::signal< void(const std::vector< std::complex< double >> &)> &eigenvalues_signal, const boost::signals2::signal< void(const FullMatrix< double > &)> &hessenberg_signal, const boost::signals2::signal< void(double)> &cond_signal)
 

Protected Attributes

AdditionalData additional_data
 
boost::signals2::signal< void(double)> condition_number_signal
 
boost::signals2::signal< void(double)> all_condition_numbers_signal
 
boost::signals2::signal< void(const std::vector< std::complex< double >> &)> eigenvalues_signal
 
boost::signals2::signal< void(const std::vector< std::complex< double >> &)> all_eigenvalues_signal
 
boost::signals2::signal< void(const FullMatrix< double > &)> hessenberg_signal
 
boost::signals2::signal< void(const FullMatrix< double > &)> all_hessenberg_signal
 
boost::signals2::signal< void(const internal::SolverGMRESImplementation::TmpVectors< VectorType > &)> krylov_space_signal
 
boost::signals2::signal< void(int)> re_orthogonalize_signal
 
FullMatrix< doubleH
 
FullMatrix< doubleH1
 
- Protected Attributes inherited from SolverBase< Vector< double > >
GrowingVectorMemory< Vector< double > > static_vector_memory
 
VectorMemory< Vector< double > > & memory
 
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const Vector< double > &current_iterate), StateCombiner > iteration_status
 

Additional Inherited Members

- Public Types inherited from SolverBase< Vector< double > >
using vector_type = Vector< double >
 

Detailed Description

template<class VectorType = Vector<double>>
class SolverGMRES< VectorType >

Implementation of the Restarted Preconditioned Direct Generalized Minimal Residual Method. The stopping criterion is the norm of the residual.

The AdditionalData structure contains the number of temporary vectors used. The size of the Arnoldi basis is this number minus three. Additionally, it allows you to choose between right or left preconditioning. The default is left preconditioning. Finally it includes a flag indicating whether or not the default residual is used as stopping criterion.

Left versus right preconditioning

AdditionalData allows you to choose between left and right preconditioning. As expected, this switches between solving for the systems P-1A and AP-1, respectively.

A second consequence is the type of residual which is used to measure convergence. With left preconditioning, this is the preconditioned residual, while with right preconditioning, it is the residual of the unpreconditioned system.

Optionally, this behavior can be overridden by using the flag AdditionalData::use_default_residual. A true value refers to the behavior described in the previous paragraph, while false reverts it. Be aware though that additional residuals have to be computed in this case, impeding the overall performance of the solver.

The size of the Arnoldi basis

The maximal basis size is controlled by AdditionalData::max_n_tmp_vectors, and it is this number minus 2. If the number of iteration steps exceeds this number, all basis vectors are discarded and the iteration starts anew from the approximation obtained so far.

Note that the minimizing property of GMRes only pertains to the Krylov space spanned by the Arnoldi basis. Therefore, restarted GMRes is not minimizing anymore. The choice of the basis length is a trade- off between memory consumption and convergence speed, since a longer basis means minimization over a larger space.

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

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.

Eigenvalue and condition number estimates

This class can estimate eigenvalues and condition number during the solution process. This is done by creating the Hessenberg matrix during the inner iterations. The eigenvalues are estimated as the eigenvalues of the Hessenberg matrix and the condition number is estimated as the ratio of the largest and smallest singular value of the Hessenberg matrix. The estimates can be obtained by connecting a function as a slot using connect_condition_number_slot and connect_eigenvalues_slot. These slots will then be called from the solver with the estimates as argument.

Author
Wolfgang Bangerth, Guido Kanschat, Ralf Hartmann.

Definition at line 178 of file solver_gmres.h.

Constructor & Destructor Documentation

◆ SolverGMRES() [1/3]

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

Constructor.

◆ SolverGMRES() [2/3]

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

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

◆ SolverGMRES() [3/3]

template<class VectorType = Vector<double>>
SolverGMRES< VectorType >::SolverGMRES ( const SolverGMRES< VectorType > &  )
delete

The copy constructor is deleted.

Member Function Documentation

◆ solve()

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

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

◆ connect_condition_number_slot()

template<class VectorType = Vector<double>>
boost::signals2::connection SolverGMRES< VectorType >::connect_condition_number_slot ( const std::function< void(double)> &  slot,
const bool  every_iteration = false 
)

Connect a slot to retrieve the estimated condition number. Called on each outer iteration if every_iteration=true, otherwise called once when iterations are ended (i.e., either because convergence has been achieved, or because divergence has been detected).

◆ connect_eigenvalues_slot()

template<class VectorType = Vector<double>>
boost::signals2::connection SolverGMRES< VectorType >::connect_eigenvalues_slot ( const std::function< void(const std::vector< std::complex< double >> &)> &  slot,
const bool  every_iteration = false 
)

Connect a slot to retrieve the estimated eigenvalues. Called on each outer iteration if every_iteration=true, otherwise called once when iterations are ended (i.e., either because convergence has been achieved, or because divergence has been detected).

◆ connect_hessenberg_slot()

template<class VectorType = Vector<double>>
boost::signals2::connection SolverGMRES< VectorType >::connect_hessenberg_slot ( const std::function< void(const FullMatrix< double > &)> &  slot,
const bool  every_iteration = true 
)

Connect a slot to retrieve the Hessenberg matrix obtained by the projection of the initial matrix on the Krylov basis. Called on each outer iteration if every_iteration=true, otherwise called once when iterations are ended (i.e., either because convergence has been achieved, or because divergence has been detected).

◆ connect_krylov_space_slot()

template<class VectorType = Vector<double>>
boost::signals2::connection SolverGMRES< VectorType >::connect_krylov_space_slot ( const std::function< void(const internal::SolverGMRESImplementation::TmpVectors< VectorType > &)> &  slot)

Connect a slot to retrieve the basis vectors of the Krylov space generated by the Arnoldi algorithm. Called at once when iterations are completed (i.e., either because convergence has been achieved, or because divergence has been detected).

◆ connect_re_orthogonalization_slot()

template<class VectorType = Vector<double>>
boost::signals2::connection SolverGMRES< VectorType >::connect_re_orthogonalization_slot ( const std::function< void(int)> &  slot)

Connect a slot to retrieve a notification when the vectors are re-orthogonalized.

◆ criterion()

template<class VectorType = Vector<double>>
virtual double SolverGMRES< VectorType >::criterion ( )
protectedvirtual

Implementation of the computation of the norm of the residual.

◆ givens_rotation()

template<class VectorType = Vector<double>>
void SolverGMRES< VectorType >::givens_rotation ( Vector< double > &  h,
Vector< double > &  b,
Vector< double > &  ci,
Vector< double > &  si,
int  col 
) const
protected

Transformation of an upper Hessenberg matrix into tridiagonal structure by givens rotation of the last column

◆ modified_gram_schmidt()

template<class VectorType = Vector<double>>
static double SolverGMRES< VectorType >::modified_gram_schmidt ( const internal::SolverGMRESImplementation::TmpVectors< VectorType > &  orthogonal_vectors,
const unsigned int  dim,
const unsigned int  accumulated_iterations,
VectorType vv,
Vector< double > &  h,
bool re_orthogonalize,
const boost::signals2::signal< void(int)> &  re_orthogonalize_signal = boost::signals2::signal< void(int)>() 
)
staticprotected

Orthogonalize the vector vv against the dim (orthogonal) vectors given by the first argument using the modified Gram-Schmidt algorithm. The factors used for orthogonalization are stored in h. The boolean re_orthogonalize specifies whether the modified Gram-Schmidt algorithm should be applied twice. The algorithm checks loss of orthogonality in the procedure every fifth step and sets the flag to true in that case. All subsequent iterations use re-orthogonalization. Calls the signal re_orthogonalize_signal if it is connected.

◆ compute_eigs_and_cond()

template<class VectorType = Vector<double>>
static void SolverGMRES< VectorType >::compute_eigs_and_cond ( const FullMatrix< double > &  H_orig,
const unsigned int  dim,
const boost::signals2::signal< void(const std::vector< std::complex< double >> &)> &  eigenvalues_signal,
const boost::signals2::signal< void(const FullMatrix< double > &)> &  hessenberg_signal,
const boost::signals2::signal< void(double)> &  cond_signal 
)
staticprotected

Estimates the eigenvalues from the Hessenberg matrix, H_orig, generated during the inner iterations. Uses these estimate to compute the condition number. Calls the signals eigenvalues_signal and cond_signal with these estimates as arguments.

Member Data Documentation

◆ additional_data

template<class VectorType = Vector<double>>
AdditionalData SolverGMRES< VectorType >::additional_data
protected

Includes the maximum number of tmp vectors.

Definition at line 320 of file solver_gmres.h.

◆ condition_number_signal

template<class VectorType = Vector<double>>
boost::signals2::signal<void(double)> SolverGMRES< VectorType >::condition_number_signal
protected

Signal used to retrieve the estimated condition number. Called once when all iterations are ended.

Definition at line 326 of file solver_gmres.h.

◆ all_condition_numbers_signal

template<class VectorType = Vector<double>>
boost::signals2::signal<void(double)> SolverGMRES< VectorType >::all_condition_numbers_signal
protected

Signal used to retrieve the estimated condition numbers. Called on each outer iteration.

Definition at line 332 of file solver_gmres.h.

◆ eigenvalues_signal

template<class VectorType = Vector<double>>
boost::signals2::signal<void(const std::vector<std::complex<double>> &)> SolverGMRES< VectorType >::eigenvalues_signal
protected

Signal used to retrieve the estimated eigenvalues. Called once when all iterations are ended.

Definition at line 339 of file solver_gmres.h.

◆ all_eigenvalues_signal

template<class VectorType = Vector<double>>
boost::signals2::signal<void(const std::vector<std::complex<double>> &)> SolverGMRES< VectorType >::all_eigenvalues_signal
protected

Signal used to retrieve the estimated eigenvalues. Called on each outer iteration.

Definition at line 346 of file solver_gmres.h.

◆ hessenberg_signal

template<class VectorType = Vector<double>>
boost::signals2::signal<void(const FullMatrix<double> &)> SolverGMRES< VectorType >::hessenberg_signal
protected

Signal used to retrieve the Hessenberg matrix. Called once when all iterations are ended.

Definition at line 352 of file solver_gmres.h.

◆ all_hessenberg_signal

template<class VectorType = Vector<double>>
boost::signals2::signal<void(const FullMatrix<double> &)> SolverGMRES< VectorType >::all_hessenberg_signal
protected

Signal used to retrieve the Hessenberg matrix. Called on each outer iteration.

Definition at line 359 of file solver_gmres.h.

◆ krylov_space_signal

template<class VectorType = Vector<double>>
boost::signals2::signal<void( const internal::SolverGMRESImplementation::TmpVectors<VectorType> &)> SolverGMRES< VectorType >::krylov_space_signal
protected

Signal used to retrieve the Krylov space basis vectors. Called once when all iterations are ended.

Definition at line 367 of file solver_gmres.h.

◆ re_orthogonalize_signal

template<class VectorType = Vector<double>>
boost::signals2::signal<void(int)> SolverGMRES< VectorType >::re_orthogonalize_signal
protected

Signal used to retrieve a notification when the vectors are re-orthogonalized.

Definition at line 373 of file solver_gmres.h.

◆ H

template<class VectorType = Vector<double>>
FullMatrix<double> SolverGMRES< VectorType >::H
protected

Projected system matrix

Definition at line 433 of file solver_gmres.h.

◆ H1

template<class VectorType = Vector<double>>
FullMatrix<double> SolverGMRES< VectorType >::H1
protected

Auxiliary matrix for inverting H

Definition at line 438 of file solver_gmres.h.


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