Dune Core Modules (2.6.0)

Dune::ArPackPlusPlus_Algorithms< BCRSMatrix, BlockVector > Class Template Reference

Wrapper to use a range of ARPACK++ eigenvalue solvers. More...

#include <dune/istl/eigenvalue/arpackpp.hh>

Public Member Functions

 ArPackPlusPlus_Algorithms (const BCRSMatrix &m, const unsigned int nIterationsMax=100000, const unsigned int verbosity_level=0)
 Construct from required parameters. More...
 
void computeSymMaxMagnitude (const Real &epsilon, BlockVector &x, Real &lambda) const
 Assume the matrix to be square, symmetric and perform IRLM to compute an approximation lambda of its dominant (i.e. largest magnitude) eigenvalue and the corresponding approximation x of an associated eigenvector. More...
 
void computeSymMinMagnitude (const Real &epsilon, BlockVector &x, Real &lambda) const
 Assume the matrix to be square, symmetric and perform IRLM to compute an approximation lambda of its least dominant (i.e. smallest magnitude) eigenvalue and the corresponding approximation x of an associated eigenvector. More...
 
void computeSymCond2 (const Real &epsilon, Real &cond_2) const
 Assume the matrix to be square, symmetric and perform IRLM to compute an approximation of its spectral condition number which, for symmetric matrices, can be expressed as the ratio of the dominant eigenvalue's magnitude and the least dominant eigenvalue's magnitude. More...
 
void computeNonSymMax (const Real &epsilon, BlockVector &x, Real &sigma) const
 Assume the matrix to be nonsymmetric and perform IRLM to compute an approximation sigma of its largest singlar value and the corresponding approximation x of an associated singular vector. More...
 
void computeNonSymMin (const Real &epsilon, BlockVector &x, Real &sigma) const
 Assume the matrix to be nonsymmetric and perform IRLM to compute an approximation sigma of its smallest singlar value and the corresponding approximation x of an associated singular vector. More...
 
void computeNonSymCond2 (const Real &epsilon, Real &cond_2) const
 Assume the matrix to be nonsymmetric and perform IRLM to compute an approximation of its spectral condition number which can be expressed as the ratio of the largest singular value and the smallest singular value. More...
 
unsigned int getIterationCount () const
 Return the number of iterations in last application of an algorithm.
 

Detailed Description

template<typename BCRSMatrix, typename BlockVector>
class Dune::ArPackPlusPlus_Algorithms< BCRSMatrix, BlockVector >

Wrapper to use a range of ARPACK++ eigenvalue solvers.

   A class template for performing some eigenvalue algorithms
   provided by the ARPACK++ library which is based on the implicitly
   restarted Arnoldi/Lanczos method (IRAM/IRLM), a synthesis of the
   Arnoldi/Lanczos process with the implicitily shifted QR technique.
   The method is designed to compute eigenvalue-eigenvector pairs of
   large scale sparse nonsymmetric/symmetric matrices. This class
   template uses the algorithms to compute the dominant (i.e. largest
   magnitude) and least dominant (i.e. smallest magnitude) eigenvalue
   as well as the spectral condition number of square, symmetric
   matrices and to compute the largest and smallest singular value as
   well as the spectral condition number of nonsymmetric matrices.
Note
For a recent version of the ARPACK++ library working with recent compiler versions see "http://reuter.mit.edu/software/arpackpatch/" or the git repository "https://github.com/m-reuter/arpackpp.git".
Note that the Arnoldi/Lanczos process currently is initialized using a vector which is randomly generated by ARPACK++. This could be changed in a future version, since ARPACK++ supports manual initialization of this vector.
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.
Todo:
Maybe make ARPACK++ parameter ncv available to the user.
Template Parameters
BCRSMatrixType of a DUNE-ISTL BCRSMatrix whose eigenvalues respectively singular values shall be considered; is assumed to have blocklevel 2.
BlockVectorType of the associated vectors; compatible with the rows of a BCRSMatrix object (if #rows >= #ncols) or its columns (if #rows < #ncols).
Author
Sebastian Westerheide.

Constructor & Destructor Documentation

◆ ArPackPlusPlus_Algorithms()

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

Construct from required parameters.

Parameters
[in]mThe DUNE-ISTL BCRSMatrix whose eigenvalues resp. singular values shall be considered.
[in]nIterationsMaxInfluences the maximum number of Arnoldi update iterations allowed; depending on the algorithm, c*nIterationsMax iterations may be performed, where c is a natural number.
[in]verbosity_levelVerbosity setting; >= 1: algorithms print a preamble and the final result, >= 2: algorithms print information about the problem solved using ARPACK++, >= 3: the final result output includes the approximated eigenvector, >= 4: sets the ARPACK(++) verbosity mode.

Member Function Documentation

◆ computeNonSymCond2()

template<typename BCRSMatrix , typename BlockVector >
void Dune::ArPackPlusPlus_Algorithms< BCRSMatrix, BlockVector >::computeNonSymCond2 ( const Real &  epsilon,
Real &  cond_2 
) const
inline

Assume the matrix to be nonsymmetric and perform IRLM to compute an approximation of its spectral condition number which can be expressed as the ratio of the largest singular value and the smallest singular value.

Parameters
[in]epsilonThe target relative accuracy of Ritz values (0 == machine precision).
[out]cond_2The approximated spectral condition number.

◆ computeNonSymMax()

template<typename BCRSMatrix , typename BlockVector >
void Dune::ArPackPlusPlus_Algorithms< BCRSMatrix, BlockVector >::computeNonSymMax ( const Real &  epsilon,
BlockVector x,
Real &  sigma 
) const
inline

Assume the matrix to be nonsymmetric and perform IRLM to compute an approximation sigma of its largest singlar value and the corresponding approximation x of an associated singular vector.

Parameters
[in]epsilonThe target relative accuracy of Ritz values (0 == machine precision).
[out]sigmaThe approximated largest singlar value.
[out]xThe associated approximated right-singular vector (if #rows >= #ncols) respectively left-singular vector (if #rows < #ncols).

◆ computeNonSymMin()

template<typename BCRSMatrix , typename BlockVector >
void Dune::ArPackPlusPlus_Algorithms< BCRSMatrix, BlockVector >::computeNonSymMin ( const Real &  epsilon,
BlockVector x,
Real &  sigma 
) const
inline

Assume the matrix to be nonsymmetric and perform IRLM to compute an approximation sigma of its smallest singlar value and the corresponding approximation x of an associated singular vector.

Parameters
[in]epsilonThe target relative accuracy of Ritz values (0 == machine precision).
[out]sigmaThe approximated smallest singlar value.
[out]xThe associated approximated right-singular vector (if #rows >= #ncols) respectively left-singular vector (if #rows < #ncols).

◆ computeSymCond2()

template<typename BCRSMatrix , typename BlockVector >
void Dune::ArPackPlusPlus_Algorithms< BCRSMatrix, BlockVector >::computeSymCond2 ( const Real &  epsilon,
Real &  cond_2 
) const
inline

Assume the matrix to be square, symmetric and perform IRLM to compute an approximation of its spectral condition number which, for symmetric matrices, can be expressed as the ratio of the dominant eigenvalue's magnitude and the least dominant eigenvalue's magnitude.

Parameters
[in]epsilonThe target relative accuracy of Ritz values (0 == machine precision).
[out]cond_2The approximated spectral condition number.

◆ computeSymMaxMagnitude()

template<typename BCRSMatrix , typename BlockVector >
void Dune::ArPackPlusPlus_Algorithms< BCRSMatrix, BlockVector >::computeSymMaxMagnitude ( const Real &  epsilon,
BlockVector x,
Real &  lambda 
) const
inline

Assume the matrix to be square, symmetric and perform IRLM to compute an approximation lambda of its dominant (i.e. largest magnitude) eigenvalue and the corresponding approximation x of an associated eigenvector.

Parameters
[in]epsilonThe target relative accuracy of Ritz values (0 == machine precision).
[out]lambdaThe approximated dominant eigenvalue.
[out]xThe associated approximated eigenvector.

◆ computeSymMinMagnitude()

template<typename BCRSMatrix , typename BlockVector >
void Dune::ArPackPlusPlus_Algorithms< BCRSMatrix, BlockVector >::computeSymMinMagnitude ( const Real &  epsilon,
BlockVector x,
Real &  lambda 
) const
inline

Assume the matrix to be square, symmetric and perform IRLM to compute an approximation lambda of its least dominant (i.e. smallest magnitude) eigenvalue and the corresponding approximation x of an associated eigenvector.

Parameters
[in]epsilonThe target relative accuracy of Ritz values (0 == machine precision).
[out]lambdaThe approximated least dominant eigenvalue.
[out]xThe associated approximated eigenvector.

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 (Nov 13, 23:29, 2024)