Reference documentation for deal.II version 9.6.0
|
#include <deal.II/lac/solver_gmres.h>
Public Member Functions | |
void | initialize (const LinearAlgebra::OrthogonalizationStrategy orthogonalization_strategy, const unsigned int max_basis_size, const bool force_reorthogonalization) |
template<typename VectorType > | |
double | orthonormalize_nth_vector (const unsigned int n, TmpVectors< VectorType > &orthogonal_vectors, const unsigned int accumulated_iterations=0, const boost::signals2::signal< void(int)> &reorthogonalize_signal=boost::signals2::signal< void(int)>()) |
const Vector< double > & | solve_projected_system (const bool orthogonalization_finished) |
const FullMatrix< double > & | get_hessenberg_matrix () const |
Private Member Functions | |
double | do_givens_rotation (const bool delayed_reorthogonalization, const int col, FullMatrix< double > &matrix, std::vector< std::pair< double, double > > &rotations, Vector< double > &rhs) |
Private Attributes | |
FullMatrix< double > | hessenberg_matrix |
FullMatrix< double > | triangular_matrix |
std::vector< std::pair< double, double > > | givens_rotations |
Vector< double > | projected_rhs |
Vector< double > | projected_solution |
Vector< double > | h |
bool | do_reorthogonalization |
LinearAlgebra::OrthogonalizationStrategy | orthogonalization_strategy |
Class that performs the Arnoldi orthogonalization process within the SolverGMRES and SolverFGMRES classes. It uses one of the algorithms in LinearAlgebra::OrthogonalizationStrategy for the work on the global vectors, transforms the resulting Hessenberg matrix into an upper triangular matrix by Givens rotations, and eventually solves the minimization problem in the projected Krylov space.
Definition at line 135 of file solver_gmres.h.
void internal::SolverGMRESImplementation::ArnoldiProcess::initialize | ( | const LinearAlgebra::OrthogonalizationStrategy | orthogonalization_strategy, |
const unsigned int | max_basis_size, | ||
const bool | force_reorthogonalization ) |
Initialize the data structures in this class with the given parameters for the solution process.
double internal::SolverGMRESImplementation::ArnoldiProcess::orthonormalize_nth_vector | ( | const unsigned int | n, |
TmpVectors< VectorType > & | orthogonal_vectors, | ||
const unsigned int | accumulated_iterations = 0, | ||
const boost::signals2::signal< void(int)> & | reorthogonalize_signal = boost::signals2::signal< void(int)>() ) |
Orthonormalize the vector at the position n
within the array orthogonal_vectors
against the n
orthonormal vectors with indices 0, ..., n - 1
using the modified or classical Gram-Schmidt algorithm. The class internally stores the factors used for orthogonalization in an upper Hessenberg matrix. For the classical Gram-Schmidt and modified Gram-Schmidt algorithms, loss of orthogonality is checked every fifth step (in case it is not yet already set via the initialize() function). In case this is detected, all subsequent iterations use re-orthogonalization as stored internally in this class, and a call to the optional signal is made.
Note that the projected Hessenberg matrix and its factorization are only consistent if n
is incremented by one for each successive call, or if n
is zero when starting to build a new orthogonal basis in restarted GMRES; an assertion will be raised if this assumption is not fulfilled.
Within this function, the factors for the QR factorization are computed alongside the Hessenberg matrix, and an estimate of the residual in the subspace is returned from this function.
const Vector< double > & internal::SolverGMRESImplementation::ArnoldiProcess::solve_projected_system | ( | const bool | orthogonalization_finished | ) |
Using the matrix and right hand side computed during the factorization, solve the underlying minimization problem for the residual in the Krylov space, returning the resulting solution as a const reference. Note that the dimension of the vector is set to the size of the Krylov space.
const FullMatrix< double > & internal::SolverGMRESImplementation::ArnoldiProcess::get_hessenberg_matrix | ( | ) | const |
Return the upper Hessenberg matrix resulting from the Gram-Schmidt orthogonalization process.
|
private |
This is a helper function to perform the incremental computation of the QR factorization of the Hessenberg matrix involved in the Arnoldi process. More precisely, it transforms the member variable hessenberg_matrix
into an upper triangular matrix R labeled matrix
, an orthogonal matrix Q represented by a vector of Givens rotations, and the associated right hand side to minimize the norm of the solution in the Krylov subspace.
More precisely, this function is called once a new column is added to the Hessenberg matrix and performs all necessary steps for that column. First, all evaluations with the Givens rotations resulting from the previous elimination steps are performed. Then, the single additional entry below the diagonal in the Hessenberg matrix is eliminated by a Givens rotation, appending a new pair of Givens factors, and the right-hand side vector in the projected system is updated. The column number col
for which the Gram-Schmidt should run needs to be given, because the algorithmic variant with delayed orthogonalization might lag by one step compared to the other sizes in the problem, and needs to perform additional computations.
In most cases, the matrices and vectors passed to this function are the member variables of the present class, but also other scenarios are supported. The function returns the modulus of the last entry in the transformed right-hand side, which is the obtained residual of the global vector x after minimization within the Krylov space.
|
private |
Projected system matrix in upper Hessenberg form.
Definition at line 200 of file solver_gmres.h.
|
private |
Upper triangular matrix that results from performing the QR factorization with Givens rotations on the upper Hessenberg matrix; the matrix Q is contained in the array givens_rotations.
Definition at line 207 of file solver_gmres.h.
|
private |
Representation of the factor Q in the QR factorization of the Hessenberg matrix.
Definition at line 213 of file solver_gmres.h.
|
private |
Right-hand side vector for orthogonalization.
Definition at line 218 of file solver_gmres.h.
|
private |
Solution vector when computing the minimization in the projected Krylov space.
Definition at line 224 of file solver_gmres.h.
|
private |
Auxiliary vector for orthogonalization.
Definition at line 229 of file solver_gmres.h.
|
private |
Flag to keep track reorthogonalization, which is checked every fifth iteration by default for LinearAlgebra::OrthogonalizationStrategy::classical_gram_schmidt and LinearAlgebra::OrthogonalizationStrategy::modified_gram_schmidt; for LinearAlgebra::OrthogonalizationStrategy::delayed_classical_gram_schmidt, no check is made.
Definition at line 239 of file solver_gmres.h.
|
private |
Selected orthogonalization algorithm.
Definition at line 244 of file solver_gmres.h.