deal.II version GIT relicensing-1936-ge55124254b 2024-10-04 21:10:00+00:00
\(\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_basis_size=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::delayed_classical_gram_schmidt)
 

Public Attributes

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

Detailed Description

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

Standardized data struct to pipe additional data to the solver.

Definition at line 363 of file solver_gmres.h.

Constructor & Destructor Documentation

◆ AdditionalData()

template<typename VectorType = Vector<double>>
SolverGMRES< VectorType >::AdditionalData::AdditionalData ( const unsigned int  max_basis_size = 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::delayed_classical_gram_schmidt 
)
explicit

Constructor. By default, set the size of the Arnoldi basis to 30. 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. Finally, the default orthogonalization algorithm is the classical Gram-Schmidt method with delayed reorthogonalization, which combines stability with fast execution, especially in parallel.

Member Data Documentation

◆ max_n_tmp_vectors

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

Maximum number of temporary vectors. Together with max_basis_size, this parameter controls the size of the Arnoldi basis, which corresponds to max_n_tmp_vectors-2 as used in previous versions of the deal.II library. SolverGMRES assumes that there are at least three temporary vectors, so this value must be greater than or equal to three. If both this variable and max_basis_size are set to a non-zero value, the choice in max_basis_size takes precedence.

Deprecated:
Use max_basis_size instead.

Definition at line 396 of file solver_gmres.h.

◆ max_basis_size

template<typename VectorType = Vector<double>>
unsigned int SolverGMRES< VectorType >::AdditionalData::max_basis_size

Maximum size of the Arnoldi basis. SolverGMRES assumes that there is at least one vector in the Arnoldi basis, so this value must be greater than or equal to one. Note that whenever this variable is set to a non-zero value, including the value set by the default constructor, this variable takes precedence over max_n_tmp_vectors.

Definition at line 405 of file solver_gmres.h.

◆ right_preconditioning

template<typename 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 414 of file solver_gmres.h.

◆ use_default_residual

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

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

Definition at line 419 of file solver_gmres.h.

◆ force_re_orthogonalization

template<typename 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 427 of file solver_gmres.h.

◆ batched_mode

template<typename 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 435 of file solver_gmres.h.

◆ orthogonalization_strategy

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

Strategy to orthogonalize vectors.

Definition at line 440 of file solver_gmres.h.


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