Dune Core Modules (unstable)
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 singular 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 singular 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
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
-
BCRSMatrix Type of a DUNE-ISTL BCRSMatrix whose eigenvalues respectively singular values shall be considered; is assumed to have blocklevel 2. BlockVector Type of the associated vectors; compatible with the rows of a BCRSMatrix object (if #rows >= #ncols) or its columns (if #rows < #ncols).
Constructor & Destructor Documentation
◆ ArPackPlusPlus_Algorithms()
|
inline |
Construct from required parameters.
- Parameters
-
[in] m The DUNE-ISTL BCRSMatrix whose eigenvalues resp. singular values shall be considered. [in] nIterationsMax Influences 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_level Verbosity 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()
|
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] epsilon The target relative accuracy of Ritz values (0 == machine precision). [out] cond_2 The approximated spectral condition number.
◆ computeNonSymMax()
|
inline |
Assume the matrix to be nonsymmetric and perform IRLM to compute an approximation sigma of its largest singular value and the corresponding approximation x of an associated singular vector.
- Parameters
-
[in] epsilon The target relative accuracy of Ritz values (0 == machine precision). [out] sigma The approximated largest singular value. [out] x The associated approximated right-singular vector (if #rows >= #ncols) respectively left-singular vector (if #rows < #ncols).
◆ computeNonSymMin()
|
inline |
Assume the matrix to be nonsymmetric and perform IRLM to compute an approximation sigma of its smallest singular value and the corresponding approximation x of an associated singular vector.
- Parameters
-
[in] epsilon The target relative accuracy of Ritz values (0 == machine precision). [out] sigma The approximated smallest singular value. [out] x The associated approximated right-singular vector (if #rows >= #ncols) respectively left-singular vector (if #rows < #ncols).
◆ computeSymCond2()
|
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] epsilon The target relative accuracy of Ritz values (0 == machine precision). [out] cond_2 The approximated spectral condition number.
◆ computeSymMaxMagnitude()
|
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] epsilon The target relative accuracy of Ritz values (0 == machine precision). [out] lambda The approximated dominant eigenvalue. [out] x The associated approximated eigenvector.
◆ computeSymMinMagnitude()
|
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] epsilon The target relative accuracy of Ritz values (0 == machine precision). [out] lambda The approximated least dominant eigenvalue. [out] x The associated approximated eigenvector.
The documentation for this class was generated from the following file:
- dune/istl/eigenvalue/arpackpp.hh
