Dune Core Modules (2.6.0)
Iterative eigenvalue algorithms based on power iteration. More...
#include <dune/istl/eigenvalue/poweriteration.hh>
Public Types | |
typedef BlockVector::field_type | Real |
Type of underlying field. | |
typedef OperatorSum | IterationOperator |
Type of iteration operator (m_ - mu_*I) | |
Public Member Functions | |
PowerIteration_Algorithms (const BCRSMatrix &m, const unsigned int nIterationsMax=1000, const unsigned int verbosity_level=0) | |
Construct from required parameters. More... | |
PowerIteration_Algorithms (const PowerIteration_Algorithms &)=delete | |
PowerIteration_Algorithms & | operator= (const PowerIteration_Algorithms &)=delete |
void | applyPowerIteration (const Real &epsilon, BlockVector &x, Real &lambda) const |
Perform the power iteration algorithm to compute an approximation lambda of the dominant (i.e. largest magnitude) eigenvalue and the corresponding approximation x of an associated eigenvector. More... | |
template<typename ISTLLinearSolver , bool avoidLinSolverCrime = false> | |
void | applyInverseIteration (const Real &epsilon, ISTLLinearSolver &solver, BlockVector &x, Real &lambda) const |
Perform the inverse iteration algorithm to compute an approximation lambda of the least dominant (i.e. smallest magnitude) eigenvalue and the corresponding approximation x of an associated eigenvector. More... | |
template<typename ISTLLinearSolver , bool avoidLinSolverCrime = false> | |
void | applyInverseIteration (const Real &gamma, const Real &epsilon, ISTLLinearSolver &solver, BlockVector &x, Real &lambda) const |
Perform the inverse iteration with shift algorithm to compute an approximation lambda of the eigenvalue closest to a given shift and the corresponding approximation x of an associated eigenvector. More... | |
template<typename ISTLLinearSolver , bool avoidLinSolverCrime = false> | |
void | applyRayleighQuotientIteration (const Real &epsilon, ISTLLinearSolver &solver, BlockVector &x, Real &lambda) const |
Perform the Rayleigh quotient iteration algorithm to compute an approximation lambda of an eigenvalue and the corresponding approximation x of an associated eigenvector. More... | |
template<typename ISTLLinearSolver , bool avoidLinSolverCrime = false> | |
void | applyTLIMEIteration (const Real &gamma, const Real &eta, const Real &epsilon, ISTLLinearSolver &solver, const Real &delta, const std::size_t &m, bool &extrnl, BlockVector &x, Real &lambda) const |
Perform the "two-level iterative method for eigenvalue calculations
(TLIME)" iteration algorithm presented in [Szyld, 1988] to compute an approximation lambda of an eigenvalue and the corresponding approximation x of an associated eigenvector. More... | |
IterationOperator & | getIterationOperator () |
Return the iteration operator (m_ - mu_*I). More... | |
const BCRSMatrix & | getIterationMatrix () const |
Return the iteration matrix (m_ - mu_*I), provided on demand when needed (e.g. for direct solvers or preconditioning). More... | |
unsigned int | getIterationCount () const |
Return the number of iterations in last application of an algorithm. | |
Protected Member Functions | |
template<typename ISTLLinearSolver > | |
void | updateShiftMu (const Real &mu, ISTLLinearSolver &solver) const |
Update shift mu_, i.e. update iteration operator/matrix (m_ - mu_*I). More... | |
Detailed Description
class Dune::PowerIteration_Algorithms< BCRSMatrix, BlockVector >
Iterative eigenvalue algorithms based on power iteration.
Given a square matrix whose eigenvalues shall be considered, this class template provides methods for performing the power iteration algorithm, the inverse iteration algorithm, the inverse iteration with shift algorithm, the Rayleigh quotient iteration algorithm and the TLIME iteration algorithm.
- Note
- Note that all algorithms except the power iteration algorithm require matrix inversion via a linear solver. When using an iterative linear solver, the algorithms become inexact "inner-outer" iterative methods. It is known that the number of inner solver iterations can increase steadily as the outer eigenvalue iteration proceeds. In this case, you should consider using a "tuned preconditioner", see e.g. [Freitag and Spence, 2008].
- In the current implementation, preconditioners like Dune::SeqILUn which are based on matrix decomposition act on the initial iteration matrix in each iteration, even for methods like the Rayleigh quotient algorithm in which the iteration matrix (m_ - mu_*I) may change in each iteration. This is due to the fact that those preconditioners currently don't support to be notified about a change of the matrix object.
- Todo:
- The current implementation is limited to DUNE-ISTL BCRSMatrix types with blocklevel 2. An extension to blocklevel >= 2 might be provided in a future version.
- Template Parameters
-
BCRSMatrix Type of a DUNE-ISTL BCRSMatrix whose eigenvalues shall be considered; is assumed to have blocklevel 2 with square blocks. BlockVector Type of the associated vectors; compatible with the rows of a BCRSMatrix object and its columns.
Constructor & Destructor Documentation
◆ PowerIteration_Algorithms() [1/2]
|
inline |
Construct from required parameters.
- Parameters
-
[in] m The square DUNE-ISTL BCRSMatrix whose eigenvalues shall be considered. [in] nIterationsMax The maximum number of iterations allowed. [in] verbosity_level Verbosity setting; >= 1: algorithms print a preamble and the final result, >= 2: algorithms print information on each iteration, >= 3: the final result output includes the approximated eigenvector.
References Dune::BCRSMatrix< B, A >::blocklevel, DUNE_THROW, Dune::BCRSMatrix< B, A >::M(), and Dune::BCRSMatrix< B, A >::N().
◆ PowerIteration_Algorithms() [2/2]
|
delete |
disallow copying (default copy constructor does a shallow copy, if copying was required a deep copy would have to be implemented due to member variables which hold a dynamically allocated object)
Member Function Documentation
◆ applyInverseIteration() [1/2]
|
inline |
Perform the inverse iteration algorithm to compute an approximation lambda of the least dominant (i.e. smallest magnitude) eigenvalue and the corresponding approximation x of an associated eigenvector.
- Template Parameters
-
ISTLLinearSolver Type of a DUNE-ISTL InverseOperator which shall be used as a linear solver. avoidLinSolverCrime The less accurate the linear solver is, the more corrupted gets the implemented computation of lambda and its associated residual. Setting this mode can help increasing their accuracy at the cost of a bit of efficiency which is beneficial e.g. when using a very inexact linear solver. Defaults to false.
- Parameters
-
[in] epsilon The target residual norm. [in] solver The DUNE-ISTL InverseOperator which shall be used as a linear solver; is assumed to be constructed using the linear operator returned by getIterationOperator() (resp. matrix returned by getIterationMatrix()). [out] lambda The approximated least dominant eigenvalue. [in,out] x The associated approximated eigenvector; shall be initialized with an estimate for an eigenvector associated with the eigenvalue which shall be approximated.
◆ applyInverseIteration() [2/2]
|
inline |
Perform the inverse iteration with shift algorithm to compute an approximation lambda of the eigenvalue closest to a given shift and the corresponding approximation x of an associated eigenvector.
- Template Parameters
-
ISTLLinearSolver Type of a DUNE-ISTL InverseOperator which shall be used as a linear solver. avoidLinSolverCrime The less accurate the linear solver is, the more corrupted gets the implemented computation of lambda and its associated residual. Setting this mode can help increasing their accuracy at the cost of a bit of efficiency which is beneficial e.g. when using a very inexact linear solver. Defaults to false.
- Parameters
-
[in] gamma The shift. [in] epsilon The target residual norm. [in] solver The DUNE-ISTL InverseOperator which shall be used as a linear solver; is assumed to be constructed using the linear operator returned by getIterationOperator() (resp. matrix returned by getIterationMatrix()). [out] lambda The approximated eigenvalue closest to gamma. [in,out] x The associated approximated eigenvector; shall be initialized with an estimate for an eigenvector associated with the eigenvalue which shall be approximated.
◆ applyPowerIteration()
|
inline |
Perform the power iteration algorithm to compute an approximation lambda of the dominant (i.e. largest magnitude) eigenvalue and the corresponding approximation x of an associated eigenvector.
- Parameters
-
[in] epsilon The target residual norm. [out] lambda The approximated dominant eigenvalue. [in,out] x The associated approximated eigenvector; shall be initialized with an estimate for an eigenvector associated with the eigenvalue which shall be approximated.
◆ applyRayleighQuotientIteration()
|
inline |
Perform the Rayleigh quotient iteration algorithm to compute an approximation lambda of an eigenvalue and the corresponding approximation x of an associated eigenvector.
- Template Parameters
-
ISTLLinearSolver Type of a DUNE-ISTL InverseOperator which shall be used as a linear solver. avoidLinSolverCrime The less accurate the linear solver is, the more corrupted gets the implemented computation of lambda and its associated residual. Setting this mode can help increasing their accuracy at the cost of a bit of efficiency which is beneficial e.g. when using a very inexact linear solver. Defaults to false.
- Parameters
-
[in] epsilon The target residual norm. [in] solver The DUNE-ISTL InverseOperator which shall be used as a linear solver; is assumed to be constructed using the linear operator returned by getIterationOperator() (resp. matrix returned by getIterationMatrix()). [in,out] lambda The approximated eigenvalue; shall be initialized with an estimate for the eigenvalue which shall be approximated. [in,out] x The associated approximated eigenvector; shall be initialized with an estimate for an eigenvector associated with the eigenvalue which shall be approximated.
◆ applyTLIMEIteration()
|
inline |
Perform the "two-level iterative method for eigenvalue calculations (TLIME)" iteration algorithm presented in [Szyld, 1988] to compute an approximation lambda of an eigenvalue and the corresponding approximation x of an associated eigenvector.
The algorithm combines the inverse iteration with shift and the Rayleigh quotient iteration in order to compute an eigenvalue in a given interval J = (gamma - eta, gamma + eta). It guarantees that if an eigenvalue exists in J, the method will converge to an eigenvalue in J, while exploiting the cubic convergence of the Rayleigh quotient iteration, but without its drawback that - depending on the initial vector - it can converge to an arbitrary eigenvalue of the matrix. When J is free of eigenvalues, the method will determine this fact and converge linearly to the eigenvalue closest to J.
- Template Parameters
-
ISTLLinearSolver Type of a DUNE-ISTL InverseOperator which shall be used as a linear solver. avoidLinSolverCrime The less accurate the linear solver is, the more corrupted gets the implemented computation of lambda and its associated residual. Setting this mode can help increasing their accuracy at the cost of a bit of efficiency which is beneficial e.g. when using a very inexact linear solver. Defaults to false.
- Parameters
-
[in] gamma An estimate for the eigenvalue which shall be approximated. [in] eta Radius around gamma in which the eigenvalue is expected. [in] epsilon The target norm of the residual with respect to the Rayleigh quotient. [in] solver The DUNE-ISTL InverseOperator which shall be used as a linear solver; is assumed to be constructed using the linear operator returned by getIterationOperator() (resp. matrix returned by getIterationMatrix()). [in] delta The target relative change of the Rayleigh quotient, indicating that inverse iteration has become stationary and switching to Rayleigh quotient iteration is appropriate; is only considered if J is free of eigenvalues. [in] m The minimum number of inverse iterations before switching to Rayleigh quotient iteration; is only considered if J is free of eigenvalues. [out] extrnl If true, the interval J is free of eigenvalues; the approximated eigenvalue-eigenvector pair (lambda,x_s) then corresponds to the eigenvalue closest to J. [out] lambda The approximated eigenvalue. [in,out] x The associated approximated eigenvector; shall be initialized with an estimate for an eigenvector associated with the eigenvalue which shall be approximated.
◆ getIterationMatrix()
|
inline |
Return the iteration matrix (m_ - mu_*I), provided on demand when needed (e.g. for direct solvers or preconditioning).
The matrix returned by this method shall be used to create the linear solver object if it requires that the matrix is provided explicitly. For linear solvers which operate completely matrix free use getIterationOperator() instead.
- Note
- Calling this method creates a new DUNE-ISTL BCRSMatrix object which requires as much memory as the matrix whose eigenvalues shall be considered.
◆ getIterationOperator()
|
inline |
Return the iteration operator (m_ - mu_*I).
The linear operator returned by this method shall be used to create the linear solver object. For linear solvers or preconditioners which require that the matrix is provided explicitly use getIterationMatrix() instead/additionally.
◆ operator=()
|
delete |
disallow copying (default assignment operator does a shallow copy, if copying was required a deep copy would have to be implemented due to member variables which hold a dynamically allocated object)
◆ updateShiftMu()
|
inlineprotected |
Update shift mu_, i.e. update iteration operator/matrix (m_ - mu_*I).
- Note
- Does nothing if new shift equals the old one.
- Template Parameters
-
ISTLLinearSolver Type of a DUNE-ISTL InverseOperator which is used as a linear solver.
- Parameters
-
[in] mu The new shift. [in] solver The DUNE-ISTL InverseOperator which is used as a linear solver.
The documentation for this class was generated from the following file:
- dune/istl/eigenvalue/poweriteration.hh