Reference documentation for deal.II version 9.3.3
|
#include <deal.II/lac/solver_gmres.h>
Classes | |
struct | AdditionalData |
Public Types | |
using | vector_type = Vector< double > |
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) |
boost::signals2::connection | connect (const std::function< SolverControl::State(const unsigned int iteration, const double check_value, const Vector< double > ¤t_iterate)> &slot) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcTooFewTmpVectors (int arg1) |
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< double > | H |
FullMatrix< double > | H1 |
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 > ¤t_iterate), StateCombiner > | iteration_status |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
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 ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
static std::mutex | mutex |
void | check_no_subscribers () const noexcept |
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.
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 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.
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.
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.
Definition at line 175 of file solver_gmres.h.
|
inherited |
SolverGMRES< VectorType >::SolverGMRES | ( | SolverControl & | cn, |
VectorMemory< VectorType > & | mem, | ||
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor.
SolverGMRES< VectorType >::SolverGMRES | ( | SolverControl & | cn, |
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor. Use an object of type GrowingVectorMemory as a default to allocate memory.
|
delete |
The copy constructor is deleted.
void SolverGMRES< VectorType >::solve | ( | const MatrixType & | A, |
VectorType & | x, | ||
const VectorType & | b, | ||
const PreconditionerType & | preconditioner | ||
) |
Solve the linear system \(Ax=b\) for x.
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).
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).
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).
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).
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.
|
protectedvirtual |
Implementation of the computation of the norm of the residual.
|
protected |
Transformation of an upper Hessenberg matrix into tridiagonal structure by givens rotation of the last column
|
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.
|
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.
|
inherited |
Connect a function object that will be called periodically within iterative solvers. This function is used to attach monitors to iterative solvers, either to determine when convergence has happened, or simply to observe the progress of an iteration. See the documentation of this class for more information.
slot | A function object specified here will, with each call, receive the number of the current iteration, the value that is used to check for convergence (typically the residual of the current iterate with respect to the linear system to be solved) and the currently best available guess for the current iterate. Note that some solvers do not update the approximate solution in every iteration but only after convergence or failure has been determined (GMRES is an example); in such cases, the vector passed as the last argument to the signal is simply the best approximate at the time the signal is called, but not the vector that will be returned if the signal's return value indicates that the iteration should be terminated. The function object must return a SolverControl::State value that indicates whether the iteration should continue, has failed, or has succeeded. The results of all connected functions will then be combined to determine what should happen with the iteration. |
|
protected |
Includes the maximum number of tmp vectors.
Definition at line 317 of file solver_gmres.h.
|
protected |
Signal used to retrieve the estimated condition number. Called once when all iterations are ended.
Definition at line 323 of file solver_gmres.h.
|
protected |
Signal used to retrieve the estimated condition numbers. Called on each outer iteration.
Definition at line 329 of file solver_gmres.h.
|
protected |
Signal used to retrieve the estimated eigenvalues. Called once when all iterations are ended.
Definition at line 336 of file solver_gmres.h.
|
protected |
Signal used to retrieve the estimated eigenvalues. Called on each outer iteration.
Definition at line 343 of file solver_gmres.h.
|
protected |
Signal used to retrieve the Hessenberg matrix. Called once when all iterations are ended.
Definition at line 349 of file solver_gmres.h.
|
protected |
Signal used to retrieve the Hessenberg matrix. Called on each outer iteration.
Definition at line 356 of file solver_gmres.h.
|
protected |
Signal used to retrieve the Krylov space basis vectors. Called once when all iterations are ended.
Definition at line 364 of file solver_gmres.h.
|
protected |
Signal used to retrieve a notification when the vectors are re-orthogonalized.
Definition at line 370 of file solver_gmres.h.
|
protected |
Projected system matrix
Definition at line 430 of file solver_gmres.h.
|
protected |
Auxiliary matrix for inverting H
Definition at line 435 of file solver_gmres.h.
|
mutableprotectedinherited |
|
protectedinherited |
|
protectedinherited |
A signal that iterative solvers can execute at the end of every iteration (or in an otherwise periodic fashion) to find out whether we should continue iterating or not. The signal may call one or more slots that each will make this determination by themselves, and the result over all slots (function calls) will be determined by the StateCombiner object.
The arguments passed to the signal are (i) the number of the current iteration; (ii) the value that is used to determine convergence (oftentimes the residual, but in other cases other quantities may be used as long as they converge to zero as the iterate approaches the solution of the linear system); and (iii) a vector that corresponds to the current best guess for the solution at the point where the signal is called. Note that some solvers do not update the approximate solution in every iteration but only after convergence or failure has been determined (GMRES is an example); in such cases, the vector passed as the last argument to the signal is simply the best approximate at the time the signal is called, but not the vector that will be returned if the signal's return value indicates that the iteration should be terminated.