Dune Core Modules (2.8.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_Algorithmsoperator= (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...
 
IterationOperatorgetIterationOperator ()
 Return the iteration operator (m_ - mu_*I). More...
 
const BCRSMatrixgetIterationMatrix () 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

template<typename BCRSMatrix, typename BlockVector>
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
BCRSMatrixType of a DUNE-ISTL BCRSMatrix whose eigenvalues shall be considered; is assumed to have blocklevel 2 with square blocks.
BlockVectorType of the associated vectors; compatible with the rows of a BCRSMatrix object and its columns.
Author
Sebastian Westerheide.

Constructor & Destructor Documentation

◆ PowerIteration_Algorithms() [1/2]

template<typename BCRSMatrix , typename BlockVector >
Dune::PowerIteration_Algorithms< BCRSMatrix, BlockVector >::PowerIteration_Algorithms ( const BCRSMatrix m,
const unsigned int  nIterationsMax = 1000,
const unsigned int  verbosity_level = 0 
)
inline

Construct from required parameters.

Parameters
[in]mThe square DUNE-ISTL BCRSMatrix whose eigenvalues shall be considered.
[in]nIterationsMaxThe maximum number of iterations allowed.
[in]verbosity_levelVerbosity 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_THROW, Dune::BCRSMatrix< B, A >::M(), and Dune::BCRSMatrix< B, A >::N().

◆ PowerIteration_Algorithms() [2/2]

template<typename BCRSMatrix , typename BlockVector >
Dune::PowerIteration_Algorithms< BCRSMatrix, BlockVector >::PowerIteration_Algorithms ( const PowerIteration_Algorithms< BCRSMatrix, BlockVector > &  )
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]

template<typename BCRSMatrix , typename BlockVector >
template<typename ISTLLinearSolver , bool avoidLinSolverCrime = false>
void Dune::PowerIteration_Algorithms< BCRSMatrix, BlockVector >::applyInverseIteration ( const Real epsilon,
ISTLLinearSolver &  solver,
BlockVector x,
Real lambda 
) const
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
ISTLLinearSolverType of a DUNE-ISTL InverseOperator which shall be used as a linear solver.
avoidLinSolverCrimeThe 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]epsilonThe target residual norm.
[in]solverThe 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]lambdaThe approximated least dominant eigenvalue.
[in,out]xThe associated approximated eigenvector; shall be initialized with an estimate for an eigenvector associated with the eigenvalue which shall be approximated.

◆ applyInverseIteration() [2/2]

template<typename BCRSMatrix , typename BlockVector >
template<typename ISTLLinearSolver , bool avoidLinSolverCrime = false>
void Dune::PowerIteration_Algorithms< BCRSMatrix, BlockVector >::applyInverseIteration ( const Real gamma,
const Real epsilon,
ISTLLinearSolver &  solver,
BlockVector x,
Real lambda 
) const
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
ISTLLinearSolverType of a DUNE-ISTL InverseOperator which shall be used as a linear solver.
avoidLinSolverCrimeThe 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]gammaThe shift.
[in]epsilonThe target residual norm.
[in]solverThe 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]lambdaThe approximated eigenvalue closest to gamma.
[in,out]xThe associated approximated eigenvector; shall be initialized with an estimate for an eigenvector associated with the eigenvalue which shall be approximated.

◆ applyPowerIteration()

template<typename BCRSMatrix , typename BlockVector >
void Dune::PowerIteration_Algorithms< BCRSMatrix, BlockVector >::applyPowerIteration ( const Real epsilon,
BlockVector x,
Real lambda 
) const
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]epsilonThe target residual norm.
[out]lambdaThe approximated dominant eigenvalue.
[in,out]xThe associated approximated eigenvector; shall be initialized with an estimate for an eigenvector associated with the eigenvalue which shall be approximated.

◆ applyRayleighQuotientIteration()

template<typename BCRSMatrix , typename BlockVector >
template<typename ISTLLinearSolver , bool avoidLinSolverCrime = false>
void Dune::PowerIteration_Algorithms< BCRSMatrix, BlockVector >::applyRayleighQuotientIteration ( const Real epsilon,
ISTLLinearSolver &  solver,
BlockVector x,
Real lambda 
) const
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
ISTLLinearSolverType of a DUNE-ISTL InverseOperator which shall be used as a linear solver.
avoidLinSolverCrimeThe 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]epsilonThe target residual norm.
[in]solverThe 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]lambdaThe approximated eigenvalue; shall be initialized with an estimate for the eigenvalue which shall be approximated.
[in,out]xThe associated approximated eigenvector; shall be initialized with an estimate for an eigenvector associated with the eigenvalue which shall be approximated.

◆ applyTLIMEIteration()

template<typename BCRSMatrix , typename BlockVector >
template<typename ISTLLinearSolver , bool avoidLinSolverCrime = false>
void Dune::PowerIteration_Algorithms< BCRSMatrix, BlockVector >::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
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
ISTLLinearSolverType of a DUNE-ISTL InverseOperator which shall be used as a linear solver.
avoidLinSolverCrimeThe 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]gammaAn estimate for the eigenvalue which shall be approximated.
[in]etaRadius around gamma in which the eigenvalue is expected.
[in]epsilonThe target norm of the residual with respect to the Rayleigh quotient.
[in]solverThe 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]deltaThe 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]mThe minimum number of inverse iterations before switching to Rayleigh quotient iteration; is only considered if J is free of eigenvalues.
[out]extrnlIf 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]lambdaThe approximated eigenvalue.
[in,out]xThe associated approximated eigenvector; shall be initialized with an estimate for an eigenvector associated with the eigenvalue which shall be approximated.

◆ getIterationMatrix()

template<typename BCRSMatrix , typename BlockVector >
const BCRSMatrix & Dune::PowerIteration_Algorithms< BCRSMatrix, BlockVector >::getIterationMatrix ( ) const
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()

template<typename BCRSMatrix , typename BlockVector >
IterationOperator & Dune::PowerIteration_Algorithms< BCRSMatrix, BlockVector >::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=()

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()

template<typename BCRSMatrix , typename BlockVector >
template<typename ISTLLinearSolver >
void Dune::PowerIteration_Algorithms< BCRSMatrix, BlockVector >::updateShiftMu ( const Real mu,
ISTLLinearSolver &  solver 
) const
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
ISTLLinearSolverType of a DUNE-ISTL InverseOperator which is used as a linear solver.
Parameters
[in]muThe new shift.
[in]solverThe DUNE-ISTL InverseOperator which is used as a linear solver.

The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 22, 23:30, 2024)