eigensolver::DirPowIt< F, G > Class Template Referenceabstract

Eigenvalue solver using the direct power iteration method. More...

#include <DirPowIt.hh>

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

Public Member Functions

virtual uint converged () const
 Returns the number of converged eigenpairs (not implemented) More...
 
virtual uint converged () const=0
 Returns the number of converged eigen pairs. More...
 
 DirPowIt (concepts::Operator< F > &A, const int nev=3, const int m=3, const int k=0, const int maxiter=300, concepts::Array< concepts::Vector< G > * > *start=0, const Real tol=0)
 Constructor. More...
 
virtual const concepts::Array< concepts::Vector< G > * > & getEF ()
 Returns an array with the eigenfunctions. More...
 
virtual const concepts::Array< G > & getEV ()
 Returns an array with the eigenvalues. More...
 
virtual uint iterations () const
 Returns the number of iterations. More...
 
virtual uint iterations () const=0
 Returns the number of iterations. More...
 
virtual ~DirPowIt ()
 Destructor. More...
 

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 sort_ (concepts::Array< Real > &)
 Sorting an array using bubble sort and returning the permutation: More...
 

Private Attributes

concepts::Operator< F > & A_
 Operator A. More...
 
bool computed_
 Are the eigenpairs already computed ? More...
 
concepts::Array< G > eigenvalues_
 Storage space for eigenvalues. More...
 
concepts::Array< concepts::Vector< G > * > eigenvectors_
 Storage space for eigenvectors. More...
 
int iter_
 Number of iterations. More...
 
int k_
 Number of highest eigenvalues to be cut off: More...
 
int m_
 Size of initial matrix V;. More...
 
int maxiter_
 Maximum number of iterations allowed. More...
 
int nev_
 Number of eigenpairs to be computed. More...
 
concepts::Array< uint > permutation_
 
concepts::Array< concepts::Vector< G > * > temp_eigenvectors_
 Temporary storage space for eigenvectors. More...
 
Real tol_
 Convergence tolerance for the eigenpairs. More...
 

Detailed Description

template<typename F, typename G>
class eigensolver::DirPowIt< F, G >

Eigenvalue solver using the direct power iteration method.

Author
Peter Kauf, 2005

Definition at line 31 of file DirPowIt.hh.

Constructor & Destructor Documentation

◆ DirPowIt()

template<typename F , typename G >
eigensolver::DirPowIt< F, G >::DirPowIt ( concepts::Operator< F > &  A,
const int  nev = 3,
const int  m = 3,
const int  k = 0,
const int  maxiter = 300,
concepts::Array< concepts::Vector< G > * > *  start = 0,
const Real  tol = 0 
)

Constructor.

Parameters
AOperator A of which we want the eigenpairs
nevNumber of eigenpairs to be computed
mNumber of columns of the iteration matrix V (must be equal to the number of columns of start)
kIt is recommended to set k = 0; k was intended to cut out the largest k from the nev eigenvalues; but to really implement such a feature one would have to do it in a more sophisticated way. Presently if k != 0 sometimes the smallest and sometimes the largest k from the nev eigenvalues are cut out. Nevertehless it seems to work for k >= 0.5*nev.
maxiterMaximum number of iterations allowed
startInitial n x m array (estimate for the eigenvectors) for iteration (n is the dimension of A and m must be bigger than nev)
tolConvergence tolerance for the eigenpairs. The default value is the machine precision for double numbers.

The routine implements the following algorithm:

n = size(A,1); V = zeros(n,m); / start
d = zeros(k,1);
for i = 1 : maxiter
  V = A*V;
  [Q,R] = qr(V,0);
  T = Q'*A*Q;
  [S,D] = eig(T);
  [l,perm] = sort(-abs(diag(D)));
  V_new = Q*S(:,perm);
  if (norm(d+l(1:k)) < tol), break; end;
  V = V_new; d = -l(1:k);
  iter = iter + 1;
end;
V = V_new(:,1:k); l = -l(1:k);

◆ ~DirPowIt()

template<typename F , typename G >
virtual eigensolver::DirPowIt< F, G >::~DirPowIt ( )
virtual

Destructor.

Member Function Documentation

◆ compute_()

template<typename F , typename G >
void eigensolver::DirPowIt< F, G >::compute_ ( )
private

Compute eigenpairs.

◆ converged() [1/2]

template<typename F , typename G >
virtual uint eigensolver::DirPowIt< F, G >::converged ( ) const
inlinevirtual

Returns the number of converged eigenpairs (not implemented)

Definition at line 85 of file DirPowIt.hh.

◆ converged() [2/2]

virtual uint eigensolver::EigenSolver< G >::converged
pure virtualinherited

Returns the number of converged eigen pairs.

◆ getEF()

template<typename F , typename G >
virtual const concepts::Array<concepts::Vector<G>*>& eigensolver::DirPowIt< F, G >::getEF ( )
virtual

Returns an array with the eigenfunctions.

Implements eigensolver::EigenSolver< G >.

◆ getEV()

template<typename F , typename G >
virtual const concepts::Array<G>& eigensolver::DirPowIt< F, G >::getEV ( )
virtual

Returns an array with the eigenvalues.

Implements eigensolver::EigenSolver< G >.

◆ info()

template<typename F , typename G >
virtual std::ostream& eigensolver::DirPowIt< F, G >::info ( std::ostream &  os) const
protectedvirtual

Returns information in an output stream.

Reimplemented from eigensolver::EigenSolver< G >.

◆ iterations() [1/2]

template<typename F , typename G >
virtual uint eigensolver::DirPowIt< F, G >::iterations ( ) const
inlinevirtual

Returns the number of iterations.

Definition at line 83 of file DirPowIt.hh.

◆ iterations() [2/2]

virtual uint eigensolver::EigenSolver< G >::iterations
pure virtualinherited

Returns the number of iterations.

◆ sort_()

template<typename F , typename G >
void eigensolver::DirPowIt< F, G >::sort_ ( concepts::Array< Real > &  )
private

Sorting an array using bubble sort and returning the permutation:

Member Data Documentation

◆ A_

template<typename F , typename G >
concepts::Operator<F>& eigensolver::DirPowIt< F, G >::A_
private

Operator A.

Definition at line 90 of file DirPowIt.hh.

◆ computed_

template<typename F , typename G >
bool eigensolver::DirPowIt< F, G >::computed_
private

Are the eigenpairs already computed ?

Definition at line 113 of file DirPowIt.hh.

◆ eigenvalues_

template<typename F , typename G >
concepts::Array<G> eigensolver::DirPowIt< F, G >::eigenvalues_
private

Storage space for eigenvalues.

Definition at line 104 of file DirPowIt.hh.

◆ eigenvectors_

template<typename F , typename G >
concepts::Array<concepts::Vector<G>*> eigensolver::DirPowIt< F, G >::eigenvectors_
private

Storage space for eigenvectors.

Definition at line 106 of file DirPowIt.hh.

◆ iter_

template<typename F , typename G >
int eigensolver::DirPowIt< F, G >::iter_
private

Number of iterations.

Definition at line 98 of file DirPowIt.hh.

◆ k_

template<typename F , typename G >
int eigensolver::DirPowIt< F, G >::k_
private

Number of highest eigenvalues to be cut off:

Definition at line 96 of file DirPowIt.hh.

◆ m_

template<typename F , typename G >
int eigensolver::DirPowIt< F, G >::m_
private

Size of initial matrix V;.

Definition at line 94 of file DirPowIt.hh.

◆ maxiter_

template<typename F , typename G >
int eigensolver::DirPowIt< F, G >::maxiter_
private

Maximum number of iterations allowed.

Definition at line 100 of file DirPowIt.hh.

◆ nev_

template<typename F , typename G >
int eigensolver::DirPowIt< F, G >::nev_
private

Number of eigenpairs to be computed.

Definition at line 92 of file DirPowIt.hh.

◆ permutation_

template<typename F , typename G >
concepts::Array<uint> eigensolver::DirPowIt< F, G >::permutation_
private

Definition at line 111 of file DirPowIt.hh.

◆ temp_eigenvectors_

template<typename F , typename G >
concepts::Array<concepts::Vector<G>*> eigensolver::DirPowIt< F, G >::temp_eigenvectors_
private

Temporary storage space for eigenvectors.

Definition at line 108 of file DirPowIt.hh.

◆ tol_

template<typename F , typename G >
Real eigensolver::DirPowIt< F, G >::tol_
private

Convergence tolerance for the eigenpairs.

Definition at line 102 of file DirPowIt.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