EigenProblemSolver Class Reference

Class to find eigenvalues and corresponding eigenvectors of a given matrix in a generalized eigenproblem, i.e. Find scalars l and non-null vectors v such that [K]{v} = l[M]{v} where [K] and [M] are symmetric matrices. The eigenproblem can be originated from a PDE. For this, we will refer to the matrices K and M as Stiffness and Mass matrices respectively. More...

#include <EigenProblemSolver.h>

Public Member Functions

 EigenProblemSolver ()
 Default constructor.
 
 EigenProblemSolver (DMatrix< real_t > &A, bool eigv=false)
 Constructor for a dense unsymmetric matrix that computes the eigenvalues.
 
 EigenProblemSolver (DSMatrix< real_t > &K, int n=0)
 Constructor for a dense symmetric matrix that computes the eigenvalues.
 
 EigenProblemSolver (SkSMatrix< real_t > &K, SkSMatrix< real_t > &M, int n=0)
 Constructor for Symmetric Skyline Matrices.
 
 EigenProblemSolver (SkSMatrix< real_t > &K, Vect< real_t > &M, int n=0)
 Constructor for Symmetric Skyline Matrices.
 
 EigenProblemSolver (DSMatrix< real_t > &A, Vect< real_t > &ev, int n=0)
 Constructor for a dense symmetric matrix that compute the eigenvalues.
 
 EigenProblemSolver (DMatrix< real_t > &A, Vect< real_t > &evr, Vect< real_t > &evi, bool eigv=false)
 Constructor for a dense matrix that compute the eigenvalues.
 
 EigenProblemSolver (Equa &eq, bool lumped=true)
 Consrtuctor using partial differential equation.
 
 ~EigenProblemSolver ()
 Destructor.
 
void set (EigenMethod m, bool sym)
 Set numerical method for eigenvalue computation.
 
void setMatrix (Matrix< real_t > *K)
 Set pointer to matrix instance.
 
void setMatrix (Matrix< real_t > *K, Matrix< real_t > *M)
 Set pointers to matrix instances.
 
void setMatrix (SkSMatrix< real_t > &K, SkSMatrix< real_t > &M)
 Set matrix instances (Symmetric matrices).
 
void setMatrix (SkSMatrix< real_t > &K, Vect< real_t > &M)
 Set matrix instances (Symmetric matrices).
 
void setMatrix (DSMatrix< real_t > &K)
 Set matrix instance (Symmetric matrix).
 
void setMatrix (DMatrix< real_t > &K)
 Set matrix instance.
 
void setNbEigv (int n)
 Set number of extracted eigenvalues.
 
void setEigenVectors ()
 Switch to compute (or not) eigenvectors as well.
 
void setPDE (Equa &eq, bool lumped=true)
 Define partial differential equation to solve.
 
int run (int nb=0, bool opt=false)
 Run the eigenproblem solver.
 
void Assembly (const Element &el, real_t *eK, real_t *eM)
 Assemble element arrays into global matrices.
 
void SAssembly (const Side &sd, real_t *sK)
 Assemble side arrays into global matrix and right-hand side.
 
int runSubSpace (size_t nb_eigv, size_t ss_dim=0)
 Run the subspace iteration solver.
 
void setSubspaceDimension (int dim)
 Define the subspace dimension.
 
void setMaxIter (int max_it)
 set maximal number of iterations.
 
void setTolerance (real_t eps)
 set tolerance value
 
int checkSturm (int &nb_found, int &nb_lost)
 Check how many eigenvalues have been found using Sturm sequence method.
 
int getNbIter () const
 Return actual number of performed iterations.
 
real_t getEigenValue (int n, int i=1) const
 Return the n-th eigenvalue.
 
void getEigenVector (int n, Vect< real_t > &v) const
 Return the n-th eigenvector.
 
void getEigenVector (int n, Vect< real_t > &vr, Vect< real_t > &vi) const
 Return the n-th eigenvector in the case of complex eigenvectors.
 

Detailed Description

Class to find eigenvalues and corresponding eigenvectors of a given matrix in a generalized eigenproblem, i.e. Find scalars l and non-null vectors v such that [K]{v} = l[M]{v} where [K] and [M] are symmetric matrices. The eigenproblem can be originated from a PDE. For this, we will refer to the matrices K and M as Stiffness and Mass matrices respectively.

Author
Rachid Touzani

Constructor & Destructor Documentation

◆ EigenProblemSolver() [1/7]

EigenProblemSolver ( DMatrix< real_t > &  A,
bool  eigv = false 
)

Constructor for a dense unsymmetric matrix that computes the eigenvalues.

This constructor solves in place the eigenvalues problem and stores them in a vector. The eigenvectors can be obtained by calling the member function getEigenVector. Eigenvalue extraction uses the reduction to Hessenberg form of the matrix.

Parameters
[in]AMatrix for which eigenmodes are sought.
[in]eigvSwitch to compute or not eigenvectors [Default: false]
Remarks
All eigenvalues are extracted.

◆ EigenProblemSolver() [2/7]

EigenProblemSolver ( DSMatrix< real_t > &  K,
int  n = 0 
)

Constructor for a dense symmetric matrix that computes the eigenvalues.

This constructor solves in place the eigenvalues problem and stores them in a vector. The eigenvectors can be obtained by calling the member function getEigenVector.

Parameters
[in]KMatrix for which eigenmodes are sought.
[in]nNumber of eigenvalues to extract. By default all eigenvalues are computed.

◆ EigenProblemSolver() [3/7]

EigenProblemSolver ( SkSMatrix< real_t > &  K,
SkSMatrix< real_t > &  M,
int  n = 0 
)

Constructor for Symmetric Skyline Matrices.

Parameters
[in]K"Stiffness" matrix
[in]M"Mass" matrix
[in]nNumber of eigenvalues to extract. By default all eigenvalues are computed.
Note
The generalized eigenvalue problem is defined by Kx = aMx, where K and M are referred to as stiffness and mass matrix.

◆ EigenProblemSolver() [4/7]

EigenProblemSolver ( SkSMatrix< real_t > &  K,
Vect< real_t > &  M,
int  n = 0 
)

Constructor for Symmetric Skyline Matrices.

Parameters
[in]K"Stiffness" matrix
[in]MDiagonal "Mass" matrix stored as a Vect instance
[in]nNumber of eigenvalues to extract. By default all eigenvalues are computed.
Note
The generalized eigenvalue problem is defined by Kx = aMx, where K and M are referred to as stiffness and mass matrix.

◆ EigenProblemSolver() [5/7]

EigenProblemSolver ( DSMatrix< real_t > &  A,
Vect< real_t > &  ev,
int  n = 0 
)

Constructor for a dense symmetric matrix that compute the eigenvalues.

This constructor solves in place the eigenvalues problem and stores them in a vector. Extraction of eigenvalues and eigenvectors uses the subspace iteration method. The eigenvectors can be obtained by calling the member function getEigenVector.

Parameters
[in]AMatrix for which eigenmodes are sought.
[in]evVector containing all computed (real) eigenvalues.
[in]nNumber of eigenvalues to extract. By default all eigenvalues are computed.
Remarks
The vector ev does not need to be sized before.

◆ EigenProblemSolver() [6/7]

EigenProblemSolver ( DMatrix< real_t > &  A,
Vect< real_t > &  evr,
Vect< real_t > &  evi,
bool  eigv = false 
)

Constructor for a dense matrix that compute the eigenvalues.

This constructor solves in place the eigenvalues problem and stores them separately in real and imaginary part vectors. Extraction of eigenvalues and eigenvectors uses the QR method after transforming the matrix in Hessenberg form. The eigenvectors can be obtained by calling the member function getEigenVector.

Parameters
[in]AMatrix for which eigenmodes are sought.
[in]evrVector containing real parts of eigenvalues.
[in]eviVector containing imaginary parts of eigenvalues.
[in]eigvSwitch to say if eigenvectors are to be computed [Default: false].
Remarks
The vectors evr and evi do not need to be sized before.

◆ EigenProblemSolver() [7/7]

EigenProblemSolver ( Equa eq,
bool  lumped = true 
)

Consrtuctor using partial differential equation.

The used equation class must have been constructed using the Mesh instance

Parameters
[in]eqReference to equation instance
[in]lumpedMass matrix is lumped (true) or not (false) [Default: true]

Member Function Documentation

◆ Assembly()

void Assembly ( const Element el,
real_t *  eK,
real_t *  eM 
)

Assemble element arrays into global matrices.

This member function is to be called from finite element equation classes

Parameters
[in]elReference to Element class
[in]eKPointer to element stiffness (or assimilated) matrix
[in]eMPointer to element mass (or assimilated) matrix

◆ checkSturm()

int checkSturm ( int &  nb_found,
int &  nb_lost 
)

Check how many eigenvalues have been found using Sturm sequence method.

Parameters
[out]nb_foundnumber of eigenvalues actually found
[out]nb_lostnumber of eigenvalues missing
Returns
- 0, Successful completion of subroutine.
  • 1, No convergent eigenvalues found.

◆ getEigenValue()

real_t getEigenValue ( int  n,
int  i = 1 
) const

Return the n-th eigenvalue.

Parameters
[in]nIndex of eigenvalue (must be > 0)
[in]ii=1: Real part, i=2: Imaginary part. [Default: 1]
Returns
Eigenvalue

◆ getEigenVector() [1/2]

void getEigenVector ( int  n,
Vect< real_t > &  v 
) const

Return the n-th eigenvector.

Parameters
[in]nLabel of eigenvector (They are stored in ascending order of eigenvalues)
[out]vVect instance where the eigenvector is stored.

◆ getEigenVector() [2/2]

void getEigenVector ( int  n,
Vect< real_t > &  vr,
Vect< real_t > &  vi 
) const

Return the n-th eigenvector in the case of complex eigenvectors.

Parameters
[in]nLabel of eigenvector (They are stored in the order of eigenvalues)
[out]vrVect instance where the real part of the eigenvector is stored.
[out]viVect instance where the imaginary part of the eigenvector is stored.

◆ run()

int run ( int  nb = 0,
bool  opt = false 
)

Run the eigenproblem solver.

Parameters
[in]nbNumber of eigenvalues to be computed. By default, all eigenvalues are computed.
[in]optToggle to compute or not corresponding eigenvectors. By default, no eigenvectors extracted.

◆ runSubSpace()

int runSubSpace ( size_t  nb_eigv,
size_t  ss_dim = 0 
)

Run the subspace iteration solver.

This function rune the Bathe subspace iteration method.

Parameters
[in]nb_eigvNumber of eigenvalues to be extracted
[in]ss_dimSubspace dimension. Must be at least equal to the number eigenvalues to seek. [Default: nb_eigv]
Returns
1: Normal execution. Convergence has been achieved. 2: Convergence for eigenvalues has not been attained.

◆ SAssembly()

void SAssembly ( const Side sd,
real_t *  sK 
)

Assemble side arrays into global matrix and right-hand side.

This member function is to be called from finite element equation classes

Parameters
[in]sdReference to Side class
[in]sKPointer to side stiffness

◆ set()

void set ( EigenMethod  m,
bool  sym 
)

Set numerical method for eigenvalue computation.

Parameters
[in]mNumerical method to compute eigenvalues:
  • SUBSPACE: Subspace iteration method (Valid for symmetric matrices only)
  • QR: QR method (Reduction to Hessenberg form)
[in]symBoolean variable to say if matrix is symmetric or not

◆ setEigenVectors()

void setEigenVectors ( )

Switch to compute (or not) eigenvectors as well.

This function is useful for the QR method. BY default eigenvectors are not extracted.

◆ setMatrix() [1/6]

void setMatrix ( DMatrix< real_t > &  K)

Set matrix instance.

This function is to be used when the default constructor is applied. Case of a standard (not generalized) eigen problem is to be solved

Parameters
[in]KMatrix instance

◆ setMatrix() [2/6]

void setMatrix ( DSMatrix< real_t > &  K)

Set matrix instance (Symmetric matrix).

This function is to be used when the default constructor is applied. Case of a standard (not generalized) eigen problem is to be solved

Parameters
[in]KMatrix instance

◆ setMatrix() [3/6]

void setMatrix ( Matrix< real_t > *  K)

Set pointer to matrix instance.

This function is to be used when the default constructor is applied. Case where the mass matrix is consistent.

Parameters
[in]KStiffness matrix pointer

◆ setMatrix() [4/6]

void setMatrix ( Matrix< real_t > *  K,
Matrix< real_t > *  M 
)

Set pointers to matrix instances.

This function is to be used when the default constructor is applied. Case where the mass matrix is consistent.

Parameters
[in]KStiffness matrix instance pointer
[in]MMass matrix instance pointer

◆ setMatrix() [5/6]

void setMatrix ( SkSMatrix< real_t > &  K,
SkSMatrix< real_t > &  M 
)

Set matrix instances (Symmetric matrices).

This function is to be used when the default constructor is applied. Case where the mass matrix is consistent.

Parameters
[in]KStiffness matrix instance
[in]MMass matrix instance

◆ setMatrix() [6/6]

void setMatrix ( SkSMatrix< real_t > &  K,
Vect< real_t > &  M 
)

Set matrix instances (Symmetric matrices).

This function is to be used when the default constructor is applied. Case where the mass matrix is (lumped) diagonal and stored in a vector.

Parameters
[in]KStiffness matrix instance
[in]MMass matrix instance where diagonal terms are stored as a vector.

◆ setNbEigv()

void setNbEigv ( int  n)

Set number of extracted eigenvalues.

This function is useful for the subspace method only.

Parameters
[in]nNumber of eigenvalues to compute. Must be at most equal to the matrix size

◆ setPDE()

void setPDE ( Equa eq,
bool  lumped = true 
)

Define partial differential equation to solve.

The used equation class must have been constructed using the Mesh instance

Parameters
[in]eqReference to equation instance
[in]lumpedMass matrix is lumped (true) or not (false) [Default: true]

◆ setSubspaceDimension()

void setSubspaceDimension ( int  dim)

Define the subspace dimension.

Parameters
[in]dimSubspace dimension. Must be larger or equal to the number of wanted eigenvalues. By default this value will be set to the number of wanted eigenvalues

◆ setTolerance()

void setTolerance ( real_t  eps)

set tolerance value

Parameters
[in]epsConvergence tolerance for eigenvalues [Default: 1.e-8]