Reference documentation for deal.II version 9.5.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\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
SolverGMRES< VectorType >::AdditionalData Struct Reference

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

Public Member Functions

 AdditionalData (const unsigned int max_n_tmp_vectors=30, const bool right_preconditioning=false, const bool use_default_residual=true, const bool force_re_orthogonalization=false, const bool batched_mode=false, const LinearAlgebra::OrthogonalizationStrategy orthogonalization_strategy=LinearAlgebra::OrthogonalizationStrategy::modified_gram_schmidt)
 

Public Attributes

unsigned int max_n_tmp_vectors
 
bool right_preconditioning
 
bool use_default_residual
 
bool force_re_orthogonalization
 
bool batched_mode
 
LinearAlgebra::OrthogonalizationStrategy orthogonalization_strategy
 

Detailed Description

template<class VectorType = Vector<double>>
struct SolverGMRES< VectorType >::AdditionalData

Standardized data struct to pipe additional data to the solver.

Definition at line 200 of file solver_gmres.h.

Constructor & Destructor Documentation

◆ AdditionalData()

template<class VectorType = Vector<double>>
SolverGMRES< VectorType >::AdditionalData::AdditionalData ( const unsigned int  max_n_tmp_vectors = 30,
const bool  right_preconditioning = false,
const bool  use_default_residual = true,
const bool  force_re_orthogonalization = false,
const bool  batched_mode = false,
const LinearAlgebra::OrthogonalizationStrategy  orthogonalization_strategy = LinearAlgebra::OrthogonalizationStrategy::modified_gram_schmidt 
)
explicit

Constructor. By default, set the number of temporary vectors to 30, i.e. do a restart every 28 iterations. Also set preconditioning from left, the residual of the stopping criterion to the default residual, and re-orthogonalization only if necessary. Also, the batched mode with reduced functionality to track information is disabled by default.

Member Data Documentation

◆ max_n_tmp_vectors

template<class VectorType = Vector<double>>
unsigned int SolverGMRES< VectorType >::AdditionalData::max_n_tmp_vectors

Maximum number of temporary vectors. This parameter controls the size of the Arnoldi basis, which for historical reasons is max_n_tmp_vectors-2. SolverGMRES assumes that there are at least three temporary vectors, so this value must be greater than or equal to three.

Definition at line 225 of file solver_gmres.h.

◆ right_preconditioning

template<class VectorType = Vector<double>>
bool SolverGMRES< VectorType >::AdditionalData::right_preconditioning

Flag for right preconditioning.

Note
Change between left and right preconditioning will also change the way residuals are evaluated. See the corresponding section in the SolverGMRES.

Definition at line 234 of file solver_gmres.h.

◆ use_default_residual

template<class VectorType = Vector<double>>
bool SolverGMRES< VectorType >::AdditionalData::use_default_residual

Flag for the default residual that is used to measure convergence.

Definition at line 239 of file solver_gmres.h.

◆ force_re_orthogonalization

template<class VectorType = Vector<double>>
bool SolverGMRES< VectorType >::AdditionalData::force_re_orthogonalization

Flag to force re-orthogonalization of orthonormal basis in every step. If set to false, the solver automatically checks for loss of orthogonality every 5 iterations and enables re-orthogonalization only if necessary.

Definition at line 247 of file solver_gmres.h.

◆ batched_mode

template<class VectorType = Vector<double>>
bool SolverGMRES< VectorType >::AdditionalData::batched_mode

Flag to control whether a reduced mode of the solver should be run. This is necessary when running (several) SolverGMRES instances involving very small and cheap linear systems where the feedback from all signals, eigenvalue computations, and log stream are disabled.

Definition at line 255 of file solver_gmres.h.

◆ orthogonalization_strategy

template<class VectorType = Vector<double>>
LinearAlgebra::OrthogonalizationStrategy SolverGMRES< VectorType >::AdditionalData::orthogonalization_strategy

Strategy to orthogonalize vectors.

Definition at line 260 of file solver_gmres.h.


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