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 Types | Public Member Functions | Public Attributes | List of all members
PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData Struct Reference

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

Public Types

enum class  EigenvalueAlgorithm { lanczos , power_iteration }
 
enum class  PolynomialType { first_kind , fourth_kind }
 

Public Member Functions

 AdditionalData (const unsigned int degree=1, const double smoothing_range=0., const unsigned int eig_cg_n_iterations=8, const double eig_cg_residual=1e-2, const double max_eigenvalue=1, const EigenvalueAlgorithm eigenvalue_algorithm=EigenvalueAlgorithm::lanczos, const PolynomialType polynomial_type=PolynomialType::first_kind)
 
AdditionalDataoperator= (const AdditionalData &other_data)
 

Public Attributes

unsigned int degree
 
double smoothing_range
 
unsigned int eig_cg_n_iterations
 
double eig_cg_residual
 
double max_eigenvalue
 
AffineConstraints< double > constraints
 
std::shared_ptr< PreconditionerType > preconditioner
 
EigenvalueAlgorithm eigenvalue_algorithm
 
PolynomialType polynomial_type
 

Detailed Description

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
struct PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData

Standardized data struct to pipe additional parameters to the preconditioner.

Definition at line 1944 of file precondition.h.

Member Enumeration Documentation

◆ EigenvalueAlgorithm

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
enum class PreconditionChebyshev::AdditionalData::EigenvalueAlgorithm
strong

An enum to define the available types of eigenvalue estimation algorithms.

Enumerator
lanczos 

This option runs the conjugate gradient solver and computes an eigenvalue estimation from the underlying Lanczos space. This only works for symmetric positive definite matrices.

power_iteration 

This option runs a power iteration to estimate the largest eigenvalue. This algorithm also works for non-symmetric matrices, but typically gives less accurate estimates than the option 'lanczos' because it does not take the relation between vectors in the iterations into account (roughly speaking the off-diagonal entries in the tri-diagonal matrix of the Lanczos iteration).

Definition at line 1950 of file precondition.h.

◆ PolynomialType

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
enum class PreconditionChebyshev::AdditionalData::PolynomialType
strong

An enum to define the available types of polynomial types.

Enumerator
first_kind 

First-kind Chebyshev polynomials.

fourth_kind 

Fourth-kind Chebyshev polynomials according to [125] and [147].

Definition at line 1972 of file precondition.h.

Constructor & Destructor Documentation

◆ AdditionalData()

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::AdditionalData ( const unsigned int  degree = 1,
const double  smoothing_range = 0.,
const unsigned int  eig_cg_n_iterations = 8,
const double  eig_cg_residual = 1e-2,
const double  max_eigenvalue = 1,
const EigenvalueAlgorithm  eigenvalue_algorithm = EigenvalueAlgorithm::lanczos,
const PolynomialType  polynomial_type = PolynomialType::first_kind 
)

Constructor.

Member Function Documentation

◆ operator=()

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
AdditionalData & PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::operator= ( const AdditionalData other_data)

Copy assignment operator.

Member Data Documentation

◆ degree

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::degree

This determines the degree of the Chebyshev polynomial. The degree of the polynomial gives the number of matrix-vector products to be performed for one application of the step() operation. During vmult(), the method performs (degree-1) matrix-vector products. Degree one corresponds to a damped Jacobi method.

If the degree is set to numbers::invalid_unsigned_int, the algorithm will automatically determine the number of necessary iterations based on the usual Chebyshev error formula as mentioned in the discussion of the main class.

Definition at line 2016 of file precondition.h.

◆ smoothing_range

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::smoothing_range

This sets the range between the largest eigenvalue in the matrix and the smallest eigenvalue to be treated. If the parameter is set to a number less than 1, an estimate for the largest and for the smallest eigenvalue will be calculated internally. For a smoothing range larger than one, the Chebyshev polynomial will act in the interval \([\lambda_\mathrm{max}/ \tt{smoothing\_range}, \lambda_\mathrm{max}]\), where \(\lambda_\mathrm{max}\) is an estimate of the maximum eigenvalue of the matrix. A choice of smoothing_range between 5 and 20 is useful in case the preconditioner is used as a smoother in multigrid.

Definition at line 2029 of file precondition.h.

◆ eig_cg_n_iterations

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::eig_cg_n_iterations

Maximum number of CG iterations performed for finding the maximum eigenvalue. If set to zero, no computations are performed. Instead, the user must supply a largest eigenvalue via the variable PreconditionChebyshev::AdditionalData::max_eigenvalue.

Definition at line 2037 of file precondition.h.

◆ eig_cg_residual

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::eig_cg_residual

Tolerance for CG iterations performed for finding the maximum eigenvalue.

Definition at line 2043 of file precondition.h.

◆ max_eigenvalue

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::max_eigenvalue

Maximum eigenvalue to work with. Only in effect if eig_cg_n_iterations is set to zero, otherwise this parameter is ignored.

Definition at line 2050 of file precondition.h.

◆ constraints

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
AffineConstraints<double> PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::constraints

Constraints to be used for the operator given. This variable is used to zero out the correct entries when creating an initial guess.

Definition at line 2056 of file precondition.h.

◆ preconditioner

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
std::shared_ptr<PreconditionerType> PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::preconditioner

Stores the preconditioner object that the Chebyshev is wrapped around.

Definition at line 2061 of file precondition.h.

◆ eigenvalue_algorithm

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
EigenvalueAlgorithm PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::eigenvalue_algorithm

Specifies the underlying eigenvalue estimation algorithm.

Definition at line 2066 of file precondition.h.

◆ polynomial_type

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
PolynomialType PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::polynomial_type

Specifies the polynomial type to be used.

Definition at line 2071 of file precondition.h.


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