Eigenvalue solver using ARPACK, the routine dnaupd or znaupd. More...

#include <ARPACK.hh>

Inheritance diagram for eigensolver::ArPack< F, G, H >:
eigensolver::EigenSolver< F > concepts::OutputOperator

Public Types

enum  modus { NORMAL = 1, REGINV = 2, SHIFTINV = 3, SHIFTINVIMAG = 4 }
 Specify mode of ARPACK which should be used to compute the Ritz values $ \nu $ of OP. More...
 
enum  which {
  LM, SM, LR, SR,
  LI, SI
}
 Specify which of the Ritz values $ \nu $ of OP (described in modus) to compute. More...
 

Public Member Functions

 ArPack (concepts::Operator< G > &OP, concepts::Operator< F > &A, concepts::Operator< H > &B, const int kmax=1, const Real tol=0.0, const int maxiter=300, enum which target=SM, enum modus mode=REGINV, const F sigma=0.0, const concepts::Vector< F > *start=0, const concepts::Array< F > *resid=0, const bool schur=false)
 Constructor. More...
 
virtual uint converged () const
 Returns the number of converged eigenpairs. More...
 
virtual const concepts::Array< concepts::Vector< F > * > & getEF ()
 Returns an array with the eigenfunctions. More...
 
virtual const concepts::Array< F > & getEV ()
 Returns an array with the real parts of the eigenvalues in the first kmax entries and the imaginary parts of the eigenvalues in the second kmax entries. More...
 
concepts::Array< F > getRESID ()
 Returns the RESID vector: More...
 
virtual uint iterations () const
 Returns the number of iterations. More...
 
virtual ~ArPack ()
 

Protected Member Functions

virtual std::ostream & info (std::ostream &os) const
 Returns information in an output stream. More...
 

Private Member Functions

void compute_ ()
 Compute eigenpairs. More...
 
void computeCmplx_ ()
 

Private Attributes

concepts::Operator< F > & A_
 Stiffness matrix. More...
 
concepts::Operator< H > & B_
 Operator B as described in modus. More...
 
bool computed_
 Eigenpairs computed? More...
 
concepts::Array< F > eigenvalues_
 Storage space for eigenvalues. More...
 
concepts::Array< F > eigenvectors_
 Storage space for eigenvectors. More...
 
concepts::Array< concepts::Vector< F > * > ev_
 References into storage for eigenvectors. More...
 
int iter_
 
int k_conv_
 
int kmax_
 Number of eigenpairs to be computed. More...
 
int maxiter_
 Maximum number of Arnoldi iterations allowed. More...
 
enum modus mode_
 Mode in which ARPACK should be used. More...
 
int numop_
 
int numopb_
 
int numreo_
 
concepts::Operator< G > & OP_
 Operator OP as described in modus. More...
 
const concepts::Array< F > * resid_in_
 Calculate vectors instead of eigenvectors Input of a residual vector (startvector) as array: More...
 
concepts::Array< F > resid_out_
 The residual vector for output: More...
 
const bool schur_
 
const F sigma_
 Shift for the shift-invert modes: More...
 
const concepts::Vector< F > * start_
 Starting vector, if not given ARPACK invents one. More...
 
enum which target_
 Sort of eigenvalues to compute. More...
 
Real tol_
 Convergence tolerance for the eigenpairs. More...
 

Detailed Description

template<typename F, typename G = Real, typename H = Real>
class eigensolver::ArPack< F, G, H >

Eigenvalue solver using ARPACK, the routine dnaupd or znaupd.

ARPACK is designed to solve large scale eigenvalue problems. The package is designed to compute a few eigenvalues and corresponding eigenvectors of a general n by n matrix A. It is most appropriate for large sparse or structured matrices A. This software is based upon an algorithmic variant of the Arnoldi process called the Implicitly Restarted Arnoldi Method (IRAM). When the matrix A is symmetric it reduces to a variant of the Lanczos process called the Implicitly Restarted Lanczos Method (IRLM). These variants may be viewed as a synthesis of the Arnoldi/Lanczos process with the Implicitly Shifted QR technique that is suitable for large scale problems.

ARPACK software is capable of solving large scale symmetric, nonsymmetric, and generalized eigenproblems from significant application areas. The software is designed to compute a few (k) eigenvalues with user specified features such as those of largest real part or largest magnitude. No auxiliary storage is required. A set of Schur basis vectors for the desired k-dimensional eigen-space is computed which is numerically orthogonal to working precision. Numerically accurate eigenvectors are available on request.

dnaupd uses implicitly restarted Arnoldi iteration to solve the generalized eigenvalue problem $ A x = \lambda B x $ with A general non-symmetric and B symmetric and positive definite. For A symmetric, use ArPackSymm instead.

znaupd uses implicitly restarted Arnoldi iteration to solve the generalized eigenvalue problem $ A x = \lambda B x $ with A general non-symmetric, complex and B symmetric and positive definite, real. For A symmetric and complex also use this routine.

See also
Richard B. Lehoucq, Kristyn J. Maschhoff, Danny C. Sorensen, and Chao Yang, ARPACK Homepage.
Richard B. Lehoucq, Danny C. Sorensen, and Chao Yang. ARPACK users' guide. Software, Environments, and Tools. Society for Industrial and Applied Mathematics (SIAM), Philadelphia, PA, 1998.
Author
Christoph Winkelmann, 2004

Definition at line 70 of file ARPACK.hh.

Member Enumeration Documentation

◆ modus

template<typename F , typename G = Real, typename H = Real>
enum eigensolver::ArPack::modus

Specify mode of ARPACK which should be used to compute the Ritz values $ \nu $ of OP.

Enumerator
NORMAL 
REGINV 
SHIFTINV 
SHIFTINVIMAG 

Definition at line 92 of file ARPACK.hh.

◆ which

template<typename F , typename G = Real, typename H = Real>
enum eigensolver::ArPack::which

Specify which of the Ritz values $ \nu $ of OP (described in modus) to compute.

Enumerator
LM 
SM 
LR 
SR 
LI 
SI 

Definition at line 75 of file ARPACK.hh.

Constructor & Destructor Documentation

◆ ArPack()

template<typename F , typename G = Real, typename H = Real>
eigensolver::ArPack< F, G, H >::ArPack ( concepts::Operator< G > &  OP,
concepts::Operator< F > &  A,
concepts::Operator< H > &  B,
const int  kmax = 1,
const Real  tol = 0.0,
const int  maxiter = 300,
enum which  target = SM,
enum modus  mode = REGINV,
const F  sigma = 0.0,
const concepts::Vector< F > *  start = 0,
const concepts::Array< F > *  resid = 0,
const bool  schur = false 
)

Constructor.

Parameters
OPOperator OP as described in modus
AStiffness matrix
BOperator B as descirbed in modus
kmaxNumber of eigenpairs to be computed
tolConvergence tolerance for the eigenpairs. The default value 0.0 is replaced by DLAMCH('EPS') from LAPACK.
maxiterMaximum number of Arnoldi iterations allowed
targetWhat sort of eigenvalues to compute
modeMode in which ARPACK should be used
sigmaShift for the shift-invert modes
startInitial vector for iteration
schurCalculate Schur vector basis instead of eigenvectors

◆ ~ArPack()

template<typename F , typename G = Real, typename H = Real>
virtual eigensolver::ArPack< F, G, H >::~ArPack ( )
virtual

Member Function Documentation

◆ compute_()

template<typename F , typename G = Real, typename H = Real>
void eigensolver::ArPack< F, G, H >::compute_ ( )
private

Compute eigenpairs.

◆ computeCmplx_()

template<typename F , typename G = Real, typename H = Real>
void eigensolver::ArPack< F, G, H >::computeCmplx_ ( )
private

◆ converged()

template<typename F , typename G = Real, typename H = Real>
virtual uint eigensolver::ArPack< F, G, H >::converged ( ) const
inlinevirtual

Returns the number of converged eigenpairs.

Implements eigensolver::EigenSolver< F >.

Definition at line 150 of file ARPACK.hh.

◆ getEF()

template<typename F , typename G = Real, typename H = Real>
virtual const concepts::Array<concepts::Vector<F>*>& eigensolver::ArPack< F, G, H >::getEF ( )
virtual

Returns an array with the eigenfunctions.

If eigenvalue k is real, entry k contains the (real) eigenfunction. If eigenvalues k and k+1 are complex conjugate, entry k contains the real part of and entry k+1 contains the imaginary part of eigenfunction k. Eigenfunction k+1 is the complex conjugate of eigenfunction k.

Implements eigensolver::EigenSolver< F >.

◆ getEV()

template<typename F , typename G = Real, typename H = Real>
virtual const concepts::Array<F>& eigensolver::ArPack< F, G, H >::getEV ( )
virtual

Returns an array with the real parts of the eigenvalues in the first kmax entries and the imaginary parts of the eigenvalues in the second kmax entries.

Implements eigensolver::EigenSolver< F >.

◆ getRESID()

template<typename F , typename G = Real, typename H = Real>
concepts::Array<F> eigensolver::ArPack< F, G, H >::getRESID ( )
inline

Returns the RESID vector:

Definition at line 152 of file ARPACK.hh.

◆ info()

template<typename F , typename G = Real, typename H = Real>
virtual std::ostream& eigensolver::ArPack< F, G, H >::info ( std::ostream &  os) const
protectedvirtual

Returns information in an output stream.

Reimplemented from eigensolver::EigenSolver< F >.

◆ iterations()

template<typename F , typename G = Real, typename H = Real>
virtual uint eigensolver::ArPack< F, G, H >::iterations ( ) const
inlinevirtual

Returns the number of iterations.

Implements eigensolver::EigenSolver< F >.

Definition at line 148 of file ARPACK.hh.

Member Data Documentation

◆ A_

template<typename F , typename G = Real, typename H = Real>
concepts::Operator<F>& eigensolver::ArPack< F, G, H >::A_
private

Stiffness matrix.

Definition at line 159 of file ARPACK.hh.

◆ B_

template<typename F , typename G = Real, typename H = Real>
concepts::Operator<H>& eigensolver::ArPack< F, G, H >::B_
private

Operator B as described in modus.

Definition at line 161 of file ARPACK.hh.

◆ computed_

template<typename F , typename G = Real, typename H = Real>
bool eigensolver::ArPack< F, G, H >::computed_
private

Eigenpairs computed?

Definition at line 188 of file ARPACK.hh.

◆ eigenvalues_

template<typename F , typename G = Real, typename H = Real>
concepts::Array<F> eigensolver::ArPack< F, G, H >::eigenvalues_
private

Storage space for eigenvalues.

Definition at line 170 of file ARPACK.hh.

◆ eigenvectors_

template<typename F , typename G = Real, typename H = Real>
concepts::Array<F> eigensolver::ArPack< F, G, H >::eigenvectors_
private

Storage space for eigenvectors.

Definition at line 172 of file ARPACK.hh.

◆ ev_

template<typename F , typename G = Real, typename H = Real>
concepts::Array<concepts::Vector<F>*> eigensolver::ArPack< F, G, H >::ev_
private

References into storage for eigenvectors.

Definition at line 174 of file ARPACK.hh.

◆ iter_

template<typename F , typename G = Real, typename H = Real>
int eigensolver::ArPack< F, G, H >::iter_
private

Definition at line 168 of file ARPACK.hh.

◆ k_conv_

template<typename F , typename G = Real, typename H = Real>
int eigensolver::ArPack< F, G, H >::k_conv_
private

Definition at line 168 of file ARPACK.hh.

◆ kmax_

template<typename F , typename G = Real, typename H = Real>
int eigensolver::ArPack< F, G, H >::kmax_
private

Number of eigenpairs to be computed.

Definition at line 163 of file ARPACK.hh.

◆ maxiter_

template<typename F , typename G = Real, typename H = Real>
int eigensolver::ArPack< F, G, H >::maxiter_
private

Maximum number of Arnoldi iterations allowed.

Definition at line 165 of file ARPACK.hh.

◆ mode_

template<typename F , typename G = Real, typename H = Real>
enum modus eigensolver::ArPack< F, G, H >::mode_
private

Mode in which ARPACK should be used.

Definition at line 174 of file ARPACK.hh.

◆ numop_

template<typename F , typename G = Real, typename H = Real>
int eigensolver::ArPack< F, G, H >::numop_
private

Definition at line 168 of file ARPACK.hh.

◆ numopb_

template<typename F , typename G = Real, typename H = Real>
int eigensolver::ArPack< F, G, H >::numopb_
private

Definition at line 168 of file ARPACK.hh.

◆ numreo_

template<typename F , typename G = Real, typename H = Real>
int eigensolver::ArPack< F, G, H >::numreo_
private

Definition at line 168 of file ARPACK.hh.

◆ OP_

template<typename F , typename G = Real, typename H = Real>
concepts::Operator<G>& eigensolver::ArPack< F, G, H >::OP_
private

Operator OP as described in modus.

Definition at line 157 of file ARPACK.hh.

◆ resid_in_

template<typename F , typename G = Real, typename H = Real>
const concepts::Array<F>* eigensolver::ArPack< F, G, H >::resid_in_
private

Calculate vectors instead of eigenvectors Input of a residual vector (startvector) as array:

Definition at line 185 of file ARPACK.hh.

◆ resid_out_

template<typename F , typename G = Real, typename H = Real>
concepts::Array<F> eigensolver::ArPack< F, G, H >::resid_out_
private

The residual vector for output:

Definition at line 190 of file ARPACK.hh.

◆ schur_

template<typename F , typename G = Real, typename H = Real>
const bool eigensolver::ArPack< F, G, H >::schur_
private

Definition at line 186 of file ARPACK.hh.

◆ sigma_

template<typename F , typename G = Real, typename H = Real>
const F eigensolver::ArPack< F, G, H >::sigma_
private

Shift for the shift-invert modes:

Definition at line 180 of file ARPACK.hh.

◆ start_

template<typename F , typename G = Real, typename H = Real>
const concepts::Vector<F>* eigensolver::ArPack< F, G, H >::start_
private

Starting vector, if not given ARPACK invents one.

Definition at line 182 of file ARPACK.hh.

◆ target_

template<typename F , typename G = Real, typename H = Real>
enum which eigensolver::ArPack< F, G, H >::target_
private

Sort of eigenvalues to compute.

Definition at line 174 of file ARPACK.hh.

◆ tol_

template<typename F , typename G = Real, typename H = Real>
Real eigensolver::ArPack< F, G, H >::tol_
private

Convergence tolerance for the eigenpairs.

Definition at line 167 of file ARPACK.hh.


The documentation for this class was generated from the following file:
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich