deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
|
#include <deal.II/lac/precondition.h>
Public Types | |
enum class | PolynomialType { first_kind , fourth_kind } |
using | EigenvalueAlgorithm = internal::EigenvalueAlgorithm |
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) | |
Public Attributes | |
unsigned int | degree |
PolynomialType | polynomial_type |
double | smoothing_range |
unsigned int | eig_cg_n_iterations |
double | eig_cg_residual |
double | max_eigenvalue |
::AffineConstraints< double > | constraints |
EigenvalueAlgorithm | eigenvalue_algorithm |
std::shared_ptr< PreconditionerType > | preconditioner |
Standardized data struct to pipe additional parameters to the preconditioner.
Definition at line 2106 of file precondition.h.
using PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::EigenvalueAlgorithm = internal::EigenvalueAlgorithm |
Definition at line 2109 of file precondition.h.
|
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 [151] and [178]. |
Definition at line 2114 of file precondition.h.
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.
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 2152 of file precondition.h.
PolynomialType PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::polynomial_type |
Specifies the polynomial type to be used.
Definition at line 2157 of file precondition.h.
|
inherited |
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 preconditioner 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 157 of file precondition.h.
|
inherited |
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 PreconditionRelaxation::AdditionalData::max_eigenvalue.
Definition at line 165 of file precondition.h.
|
inherited |
Tolerance for iterations performed for finding the maximum eigenvalue by the eigenvalue algorithm (Lanczos or power iteration).
Definition at line 171 of file precondition.h.
|
inherited |
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 178 of file precondition.h.
|
inherited |
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 184 of file precondition.h.
|
inherited |
Stores the preconditioner object that the Chebyshev is wrapped around.
Definition at line 189 of file precondition.h.
|
inherited |
Preconditioner.
Definition at line 194 of file precondition.h.