Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/solver_qmrs.h>
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 &preconditioner) |
virtual void | print_vectors (const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const |
Public Member Functions inherited from SolverBase< VectorType > | |
SolverBase (SolverControl &solver_control, VectorMemory< VectorType > &vector_memory) | |
SolverBase (SolverControl &solver_control) | |
boost::signals2::connection | connect (const std::function< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType ¤t_iterate)> &slot) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (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) |
Protected Attributes | |
AdditionalData | additional_data |
Protected Attributes inherited from SolverBase< VectorType > | |
GrowingVectorMemory< VectorType > | static_vector_memory |
VectorMemory< VectorType > & | memory |
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType ¤t_iterate), StateCombiner > | iteration_status |
Private Member Functions | |
template<typename MatrixType , typename PreconditionerType > | |
IterationResult | iterate (const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner, VectorType &r, VectorType &u, VectorType &q, VectorType &t, VectorType &d) |
Private Attributes | |
unsigned int | step |
Additional Inherited Members | |
Public Types inherited from SolverBase< VectorType > | |
using | vector_type = VectorType |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
The SQMR (symmetric quasi-minimal residual) method is supposed to solve symmetric indefinite linear systems with symmetric, not necessarily definite preconditioners. It is a variant of the original quasi-minimal residual method (QMR) and produces the same iterative solution. This version of SQMR is adapted from the respective symmetric QMR-from-BiCG algorithm given by both Freund/Nachtigal: A new Krylov-subspace method for symmetric indefinite linear systems, NASA STI/Recon Technical Report N, 95 (1994) and Freund/Nachtigal: Software for simplified Lanczos and QMR algorithms, Appl. Num. Math. 19 (1995), pp. 319-341 and provides both right and left (but not split) preconditioning.
Note, that the QMR implementation that the given algorithm is based on is derived from classical BiCG. It can be shown (Freund/Szeto: A transpose-free quasi-minimal residual squared algorithm for non-Hermitian linear systems, Advances in Computer Methods for Partial Differential Equations VII (IMACS, New Brunswick, NJ, 1992) pp. 258-264) that the QMR iterates can be generated from the BiCG iteration through one additional vector and some scalar updates. Possible breakdowns (or precisely, divisions by zero) of BiCG therefore obviously transfer to this simple no-look-ahead algorithm.
In return the algorithm is cheap compared to classical QMR or BiCGStab, using only one matrix-vector product with the system matrix and one application of the preconditioner per iteration respectively.
The residual used for measuring convergence is only approximately calculated by an upper bound. If this value comes below a threshold prescribed within the AdditionalData struct, then the exact residual of the current QMR iterate will be calculated using another multiplication with the system matrix. By experience (according to Freund and Nachtigal) this technique is useful for a threshold that is ten times the solving tolerance, and in that case will be only used in the last one or two steps of the complete iteration.
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 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 95 of file solver_qmrs.h.
SolverQMRS< VectorType >::SolverQMRS | ( | SolverControl & | cn, |
VectorMemory< VectorType > & | mem, | ||
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor.
SolverQMRS< VectorType >::SolverQMRS | ( | SolverControl & | cn, |
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor. Use an object of type GrowingVectorMemory as a default to allocate memory.
void SolverQMRS< VectorType >::solve | ( | const MatrixType & | A, |
VectorType & | x, | ||
const VectorType & | b, | ||
const PreconditionerType & | preconditioner | ||
) |
Solve the linear system \(Ax=b\) for x.
|
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.
|
private |
The iteration loop itself. The function returns a structure indicating what happened in this function.
|
protected |
Additional parameters.
Definition at line 200 of file solver_qmrs.h.
|
private |
Number of the current iteration (accumulated over restarts)
Definition at line 235 of file solver_qmrs.h.