Public Types | Public Member Functions | Protected Types | Protected Attributes | Friends | List of all members
Seldon::EigenProblem_Base< T, MatStiff, MatMass > Class Template Referenceabstract

Base class to solve an eigenvalue problem. More...

#include <EigenvalueSolver.hxx>

Public Types

enum  {
  REGULAR_MODE, SHIFTED_MODE, IMAG_SHIFTED_MODE, INVERT_MODE,
  BUCKLING_MODE, CAYLEY_MODE
}
 several available modes to find eigenvalues (Arpack) More...
 
enum  { SMALL_EIGENVALUES, LARGE_EIGENVALUES, CENTERED_EIGENVALUES }
 parts of the spectrum (near from 0, at infinity or around a given value) More...
 
enum  { SORTED_REAL, SORTED_IMAG, SORTED_MODULUS }
 different sorting strategies
 
enum  { SOLVER_LOBPCG, SOLVER_BKS, SOLVER_BD }
 different solvers (Anasazi) More...
 
enum  { ORTHO_DGKS, ORTHO_SVQB }
 orthogonalization managers (Anasazi)
 
enum  {
  REGULAR_MODE, SHIFTED_MODE, IMAG_SHIFTED_MODE, INVERT_MODE,
  BUCKLING_MODE, CAYLEY_MODE
}
 several available modes to find eigenvalues (Arpack) More...
 
enum  { REAL_PART, IMAG_PART, COMPLEX_PART }
 real part, imaginary part or complex solution
 
typedef MatMass::value_type MassValue
 type for number stored in mass matrix
 

Public Member Functions

 EigenProblem_Base ()
 default constructor
 
void Init (int n)
 initialisation of the size of the eigenvalue problem
 
void InitMatrix (MatStiff &)
 initialization of a standard eigenvalue problem More...
 
void InitMatrix (MatStiff &, MatMass &)
 initialization of a generalized eigenvalue problem More...
 
int GetComputationalMode () const
 returns the spectral transformation used for evaluation of eigenvalues
 
void SetComputationalMode (int mode)
 sets the spectral transformation used for evaluation of eigenvalues
 
int GetNbAskedEigenvalues () const
 returns the number of eigenvalues asked by the user
 
void SetNbAskedEigenvalues (int n)
 sets the number of eigenvalues to compute
 
int GetNbAdditionalEigenvalues () const
 returns the additional number of eigenvalues
 
void SetNbAdditionalEigenvalues (int n)
 sets the number of additional eigenvalues
 
int GetNbBlocks () const
 returns the number of blocks used in blocked solvers
 
void SetNbBlocks (int)
 returns the number of blocks used in blocked solvers
 
int GetNbMaximumRestarts () const
 returns the restart parameter used in blocked solvers
 
void SetNbMaximumRestarts (int)
 sets the restart parameter used in blocked solvers
 
int GetOrthoManager () const
 returns orthogonalization manager set in Anasazi
 
int GetEigensolverType () const
 returns the solver used in Anasazi
 
void SetEigensolverType (int type)
 sets the solver used in Anasazi
 
int GetTypeSpectrum () const
 returns the spectrum desired (large, small eigenvalues, etc)
 
int GetTypeSorting () const
 returns how eigenvalues are sorted (real, imaginary part or modulus)
 
GetShiftValue () const
 returns the shift value used More...
 
GetImagShiftValue () const
 returns the imaginary part of shift value used More...
 
void SetShiftValue (const T &)
 Sets the real part of shift value.
 
void SetImagShiftValue (const T &)
 Sets the imaginary part of shift value.
 
void SetTypeSpectrum (int type, const T &val, int type_sort=SORTED_MODULUS)
 sets which eigenvalues are searched More...
 
void SetTypeSpectrum (int type, const complex< T > &val, int type_sort=SORTED_MODULUS)
 sets which eigenvalues are searched More...
 
double GetLowerBoundInterval () const
 returns lower bound of the interval where eigenvalues are searched
 
double GetUpperBoundInterval () const
 returns upper bound of the interval where eigenvalues are searched
 
void SetIntervalSpectrum (double, double)
 sets the interval where eigenvalues are searched
 
void SetCholeskyFactoForMass (bool chol=true)
 indicates the use of Cholesky factorisation in order to solve a standard eigenvalue problem instead of a generalized one
 
bool UseCholeskyFactoForMass () const
 returns true if Cholesky factorisation has to be used for mass matrix
 
void SetDiagonalMass (bool diag=true)
 indicates that the mass matrix is diagonal
 
bool DiagonalMass () const
 returns true if the mass matrix is diagonal
 
void SetStoppingCriterion (double eps)
 modifies the stopping critertion
 
double GetStoppingCriterion () const
 returns the stopping criterion
 
void SetNbMaximumIterations (int n)
 sets the maximal number of iterations allowed for the iterative algorithm
 
int GetNbMaximumIterations () const
 returns the maximal number of iterations allowed for the iterative algorithm
 
int GetNbMatrixVectorProducts () const
 returns the number of matrix-vector products performed since last call to Init
 
int GetNbArnoldiVectors () const
 returns the number of Arnoldi vectors to use
 
void SetNbArnoldiVectors (int n)
 sets the number of Arnoldi vectors to use
 
int GetM () const
 returns number of rows
 
int GetN () const
 returns number of columns
 
int GetPrintLevel () const
 returns level of verbosity
 
void SetPrintLevel (int lvl)
 sets the level of verbosity
 
void IncrementProdMatVect ()
 increment of the number of matrix vector products
 
void PrintErrorInit () const
 prints error of initialization and aborts program
 
bool IsSymmetricProblem () const
 returns true if the matrix is symmetric
 
bool IsHermitianProblem () const
 returns true if the matrix is hermitian
 
void FactorizeDiagonalMass ()
 computation of D^1/2 from D
 
template<class T0 >
void MltInvSqrtDiagonalMass (Vector< T0 > &X)
 multiplication of X by D^-1/2
 
template<class T0 >
void MltSqrtDiagonalMass (Vector< T0 > &X)
 multiplication of X by D^1/2
 
void ComputeDiagonalMass ()
 computation of diagonal of mass matrix
 
void ComputeMassForCholesky ()
 computation of mass matrix
 
void ComputeMassMatrix ()
 computation of mass matrix M
 
void MltMass (const Vector< T > &X, Vector< T > &Y)
 matrix vector product with mass matrix Y = M X
 
void ComputeStiffnessMatrix ()
 computation of stiffness matrix K
 
void ComputeStiffnessMatrix (const T &a, const T &b)
 computation of matrix a M + b*K
 
void MltStiffness (const Vector< T > &X, Vector< T > &Y)
 matrix vector product with stifness matrix Y = K X
 
void MltStiffness (const T &a, const T &b, const Vector< T > &X, Vector< T > &Y)
 matrix vector product with stifness and mass matrix Y = (a M + b K) X
 
void ComputeAndFactorizeStiffnessMatrix (const T &a, const T &b)
 computation of matrix a M + b K and factorisation of this matrix More...
 
void ComputeAndFactorizeStiffnessMatrix (const complex< T > &a, const complex< T > &b, bool real_p=true)
 computation of matrix a M + b K and factorisation of this matrix More...
 
template<class T0 >
void ComputeSolution (const Vector< T0 > &X, Vector< T0 > &Y)
 solving the linear system (a M + b K) Y = X
 
template<class TransA , class T0 >
void ComputeSolution (const TransA &transA, const Vector< T0 > &X, Vector< T0 > &Y)
 solving the linear system (a M + b K) Y = X
 
void FactorizeCholeskyMass ()
 computation of Cholesky factorisation of M from matrix M
 
template<class TransStatus >
void MltCholeskyMass (const TransStatus &transA, Vector< T > &X)
 computation of L X or L^T x if M = L L^T
 
template<class TransStatus >
void SolveCholeskyMass (const TransStatus &transA, Vector< T > &X)
 computation of L^-1 X or L^-T x if M = L L^T
 
void Clear ()
 memory release
 
int GetComputationalMode () const
 
void SetComputationalMode (int mode)
 
int GetNbAdditionalEigenvalues () const
 
void SetNbAdditionalEigenvalues (int n)
 
AnasaziParamGetAnasaziParameters ()
 returns parameters specific to Anasazi
 
SlepcParamGetSlepcParameters ()
 returns parameters specific to Slepc
 
FeastParamGetFeastParameters ()
 returns parameters specific to Feast
 
void SetCholeskyFactoForMass (bool chol=true)
 
bool UseCholeskyFactoForMass () const
 
void SetDiagonalMass (bool diag=true)
 
bool DiagonalMass () const
 
void PrintErrorInit () const
 
virtual void ComputeDiagonalMass ()=0
 
virtual void FactorizeDiagonalMass ()=0
 
virtual void GetSqrtDiagonal (Vector< T > &)=0
 
virtual void MltInvSqrtDiagonalMass (Vector< Treal > &X)=0
 
virtual void MltSqrtDiagonalMass (Vector< Treal > &X)=0
 
virtual void MltInvSqrtDiagonalMass (Vector< Tcplx > &X)=0
 
virtual void MltSqrtDiagonalMass (Vector< Tcplx > &X)=0
 
virtual void ComputeMassForCholesky ()
 
virtual void ComputeMassMatrix ()
 
virtual void MltMass (const Vector< Treal > &X, Vector< Treal > &Y)=0
 
virtual void MltMass (const Vector< Tcplx > &X, Vector< Tcplx > &Y)=0
 
virtual void MltMass (const SeldonTranspose &, const Vector< Treal > &X, Vector< Treal > &Y)=0
 
virtual void MltMass (const SeldonTranspose &, const Vector< Tcplx > &X, Vector< Tcplx > &Y)=0
 
virtual void ComputeStiffnessMatrix ()
 
virtual void ComputeStiffnessMatrix (const T &a, const T &b)
 
virtual void MltStiffness (const Vector< Treal > &X, Vector< Treal > &Y)=0
 
virtual void MltStiffness (const Vector< Tcplx > &X, Vector< Tcplx > &Y)=0
 
virtual void MltStiffness (const T &a, const T &b, const Vector< Treal > &X, Vector< Treal > &Y)=0
 
virtual void MltStiffness (const T &a, const T &b, const Vector< Tcplx > &X, Vector< Tcplx > &Y)=0
 
virtual void MltStiffness (const SeldonTranspose &, const Vector< Treal > &X, Vector< Treal > &Y)=0
 
virtual void MltStiffness (const SeldonTranspose &, const Vector< Tcplx > &X, Vector< Tcplx > &Y)=0
 
virtual void ComputeAndFactorizeStiffnessMatrix (const Treal &a, const Treal &b, int real_p=COMPLEX_PART)
 computation of matrix a M + b K and factorisation of this matrix More...
 
virtual void ComputeAndFactorizeStiffnessMatrix (const Tcplx &a, const Tcplx &b, int real_p=COMPLEX_PART)
 computation of matrix a M + b K and factorisation of this matrix More...
 
virtual void ComputeSolution (const Vector< Treal > &X, Vector< Treal > &Y)
 solving the linear system (a M + b K) Y = X
 
virtual void ComputeSolution (const Vector< Tcplx > &X, Vector< Tcplx > &Y)
 solving the linear system (a M + b K) Y = X
 
virtual void ComputeSolution (const SeldonTranspose &, const Vector< Treal > &X, Vector< Treal > &Y)
 solving the linear system (a M + b K) Y = X
 
virtual void ComputeSolution (const SeldonTranspose &, const Vector< Tcplx > &X, Vector< Tcplx > &Y)
 solving the linear system (a M + b K) Y = X
 
virtual void FactorizeCholeskyMass ()
 
virtual void MltCholeskyMass (const SeldonTranspose &transA, Vector< Treal > &X)
 computation of L X or L^T x if M = L L^T
 
virtual void MltCholeskyMass (const SeldonTranspose &transA, Vector< Tcplx > &X)
 computation of L X or L^T x if M = L L^T
 
virtual void SolveCholeskyMass (const SeldonTranspose &transA, Vector< Treal > &X)
 computation of L^-1 X or L^-T x if M = L L^T
 
virtual void SolveCholeskyMass (const SeldonTranspose &transA, Vector< Tcplx > &X)
 computation of L^-1 X or L^-T x if M = L L^T
 
virtual void Clear ()
 

Protected Types

typedef ClassComplexType< T >::Tcplx Tcplx
 
typedef ClassComplexType< T >::Treal Treal
 

Protected Attributes

int eigenvalue_computation_mode
 mode used to find eigenvalues (regular, shifted, Cayley, etc)
 
int nb_eigenvalues_wanted
 number of eigenvalues to be computed
 
int nb_add_eigenvalues
 additional number of eigenvalues More...
 
int type_spectrum_wanted
 which spectrum ? Near from Zero ? Near from Infinity ? or near from a value ?
 
int type_sort_eigenvalues
 large eigenvalues because of their real part, imaginary part or magnitude ?
 
bool use_cholesky
 if true, the generalized problem is reduced to a standard problem More...
 
bool diagonal_mass
 if true, the generalized problem is reduced to a standard problem More...
 
double stopping_criterion
 threshold for Arpack's iterative process
 
int nb_maximum_iterations
 Maximal number of iterations.
 
int nb_prod
 number of matrix-vector products
 
int n_
 size of the problem
 
shift
 shift sigma (if type_spectrum = centered_eigenvalues)
 
shift_imag
 
int nb_arnoldi_vectors
 number of Arnoldi vectors
 
bool automatic_selection_arnoldi_vectors
 if true nb_arnoldi_vectors is automatically computed
 
int print_level
 print level
 
Vector< MassValuesqrt_diagonal_mass
 diagonal D^1/2 if the mass matrix is diagonal positive
 
bool complex_system
 if true consider Real( (a M + bK)^-1) or Imag( (a M + b K)^-1 ) or the whole system, a and/or b being complex
 
MatMass * Mh
 mass matrix
 
MatStiff * Kh
 stiffness matrix
 
int type_solver
 which solver ?
 
int ortho_manager
 orthogonalization manager
 
int nb_blocks
 number of blocks for blocked solvers
 
int restart_number
 restart parameter for blocked solvers
 
double emin_interval
 interval where eigenvalues are searched
 
double emax_interval
 
int selected_part
 
AnasaziParam anasazi_param
 additional parameters for Anasazi
 
SlepcParam slepc_param
 additional parameters for Slepc
 
FeastParam feast_param
 additional parameters for Feast
 

Friends

template<class EigenPb , class Vector1 , class Matrix1 , class T0 >
void ApplyScalingEigenvec (EigenPb &var, Vector1 &eigen_values, Vector1 &lambda_imag, Matrix1 &eigen_vectors, const T0 &shiftr, const T0 &shifti)
 modification of eigenvectors to take into account scaling by mass matrix More...
 
template<class EigenPb , class Vector1 , class Matrix1 , class T0 >
void ApplyScalingEigenvec (EigenPb &var, Vector1 &eigen_values, Vector1 &lambda_imag, Matrix1 &eigen_vectors, const complex< T0 > &shiftr, const complex< T0 > &shifti)
 modification of eigenvectors to take into account scaling by mass matrix More...
 
template<class T0 , class Prop , class Storage >
void ApplyScalingEigenvec (EigenProblem_Base< T0 > &var, Vector< T0 > &eigen_values, Vector< T0 > &lambda_imag, Matrix< T0, Prop, Storage > &eigen_vectors, const T0 &shiftr, const T0 &shifti)
 modification of eigenvectors to take into account scaling by mass matrix More...
 
template<class T0 , class Prop , class Storage >
void ApplyScalingEigenvec (EigenProblem_Base< complex< T0 > > &var, Vector< complex< T0 > > &eigen_values, Vector< complex< T0 > > &lambda_imag, Matrix< complex< T0 >, Prop, Storage > &eigen_vectors, const complex< T0 > &shiftr, const complex< T0 > &shifti)
 modification of eigenvectors to take into account scaling by mass matrix More...
 

Detailed Description

template<class T, class MatStiff = Matrix<T>, class MatMass = Matrix<double>>
class Seldon::EigenProblem_Base< T, MatStiff, MatMass >

Base class to solve an eigenvalue problem.

Base class to solve a linear eigenvalue problem.

Resolution of a standard eigenvalue problem : K x = lambda x or a generalized eigenvalue problem K x = lambda M x M is called mass matrix, K stiffness matrix, lambda is the eigenvalue and x the eigenvector.

This class should not be instantiated directly, but rather derived classes like DenseEigenProblem, SparseEigenProblem, MatrixFreeEigenProblem.

Solution of a standard eigenvalue problem : K x = lambda x or a generalized eigenvalue problem K x = lambda M x M is called mass matrix, K stiffness matrix, lambda is the eigenvalue and x the eigenvector.

This class should not be instantiated directly, but rather derived classes like DenseEigenProblem, SparseEigenProblem, VirtualEigenProblem.

Definition at line 18 of file EigenvalueSolver.hxx.

Member Enumeration Documentation

◆ anonymous enum

template<class T , class MatStiff = Matrix<T>, class MatMass = Matrix<double>>
anonymous enum

several available modes to find eigenvalues (Arpack)

REGULAR_MODE : Regular mode standard problem => no linear system to solve generalized problem => M^-1 K x = lambda x (inverse of M required) SHIFTED_MODE : Shifted mode standard problem => (K - sigma I)^-1 X = lambda X generalized problem => (K - sigma M)^-1 M X = lambda X BUCKLING_MODE : Buckling mode (real symmetric problem) generalized problem => (K - sigma M)^-1 K X = lambda X CAYLEY_MODE : Cayley mode (real symmetric problem) generalized problem => (K - sigma M)^-1 (K + sigma M) X = lambda X INVERT_MODE : Shifted mode on matrix M^-1 K
IMAG_SHIFTED_MODE : using Imag( (K - sigma M)^-1 M) instead of Real( (K - sigma M)^-1 M) mode 4 in Arpack (dnaupd.f)

Definition at line 441 of file VirtualEigenvalueSolver.hxx.

◆ anonymous enum

template<class T , class MatStiff = Matrix<T>, class MatMass = Matrix<double>>
anonymous enum

several available modes to find eigenvalues (Arpack)

REGULAR_MODE : Regular mode standard problem => no linear system to solve generalized problem => M^-1 K x = lambda x (inverse of M required) SHIFTED_MODE : Shifted mode standard problem => (K - sigma I)^-1 X = lambda X generalized problem => (K - sigma M)^-1 M X = lambda X BUCKLING_MODE : Buckling mode (real symmetric problem) generalized problem => (K - sigma M)^-1 K X = lambda X CAYLEY_MODE : Cayley mode (real symmetric problem) generalized problem => (K - sigma M)^-1 (K + sigma M) X = lambda X INVERT_MODE : Shifted mode on matrix M^-1 K
IMAG_SHIFTED_MODE : using Imag( (K - sigma M)^-1 M) instead of Real( (K - sigma M)^-1 M) mode 4 in Arpack (dnaupd.f)

Definition at line 38 of file EigenvalueSolver.hxx.

◆ anonymous enum

template<class T , class MatStiff = Matrix<T>, class MatMass = Matrix<double>>
anonymous enum

parts of the spectrum (near from 0, at infinity or around a given value)

SMALL_EIGENVALUES : seeking eigenvalues near 0 LARGE_EIGENVALUES : seeking largest eigenvalues CENTERED_EIGENVALUES : seeking eigenvalues near the shift sigma

Definition at line 47 of file EigenvalueSolver.hxx.

◆ anonymous enum

template<class T , class MatStiff = Matrix<T>, class MatMass = Matrix<double>>
anonymous enum

different solvers (Anasazi)

SOLVER_LOBPCG : Locally Optimal Block Preconditioned Conjugate Gradient SOLVER_BKS : Block Krylov Schur SOLVER_BD : Block Davidson

Definition at line 58 of file EigenvalueSolver.hxx.

Member Function Documentation

◆ ComputeAndFactorizeStiffnessMatrix() [1/4]

template<class T , class MatStiff , class MatMass >
void Seldon::EigenProblem_Base< T, MatStiff, MatMass >::ComputeAndFactorizeStiffnessMatrix ( const complex< T > &  a,
const complex< T > &  b,
bool  real_part = true 
)

computation of matrix a M + b K and factorisation of this matrix

The factorisation process can be also the construction of preconditioning if an iterative solver is used to solve linear system (a M + b K) y = x

Definition at line 771 of file EigenvalueSolver.cxx.

◆ ComputeAndFactorizeStiffnessMatrix() [2/4]

template<class T , class MatStiff , class MatMass >
void Seldon::EigenProblem_Base< T, MatStiff, MatMass >::ComputeAndFactorizeStiffnessMatrix ( const T &  a,
const T &  b 
)

computation of matrix a M + b K and factorisation of this matrix

The factorisation process can be also the construction of preconditioning if an iterative solver is used to solve linear system (a M + b K) y = x

Definition at line 758 of file EigenvalueSolver.cxx.

◆ ComputeAndFactorizeStiffnessMatrix() [3/4]

template<class T >
void Seldon::EigenProblem_Base< T >::ComputeAndFactorizeStiffnessMatrix ( const Tcplx &  a,
const Tcplx &  b,
int  which_part = COMPLEX_PART 
)
virtual

computation of matrix a M + b K and factorisation of this matrix

The factorisation process can be also the construction of preconditioning if an iterative solver is used to solve linear system (a M + b K) y = x

Reimplemented in Seldon::DenseEigenProblem< T, Prop, Storage, Tmass, PropM, StorageM >, and Seldon::SparseEigenProblem< T, MatStiff, MatMass >.

Definition at line 941 of file VirtualEigenvalueSolver.cxx.

◆ ComputeAndFactorizeStiffnessMatrix() [4/4]

template<class T >
void Seldon::EigenProblem_Base< T >::ComputeAndFactorizeStiffnessMatrix ( const Treal &  a,
const Treal &  b,
int  which_part = COMPLEX_PART 
)
virtual

computation of matrix a M + b K and factorisation of this matrix

The factorisation process can be also the construction of preconditioning if an iterative solver is used to solve linear system (a M + b K) y = x

Reimplemented in Seldon::DenseEigenProblem< T, Prop, Storage, Tmass, PropM, StorageM >, and Seldon::SparseEigenProblem< T, MatStiff, MatMass >.

Definition at line 927 of file VirtualEigenvalueSolver.cxx.

◆ GetImagShiftValue()

template<class T , class MatStiff , class MatMass >
T Seldon::EigenProblem_Base< T, MatStiff, MatMass >::GetImagShiftValue

returns the imaginary part of shift value used

If type_spectrum_wanted is set to CENTERED_EIGENVALUES, we search closest eigenvalues to the shift value. Matrix (A - (shift + i shift_imag)*I)^{-1} will be used instead of A shift_imag is accessed only for real unsymmetric problems

Definition at line 288 of file EigenvalueSolver.cxx.

◆ GetShiftValue()

template<class T , class MatStiff , class MatMass >
T Seldon::EigenProblem_Base< T, MatStiff, MatMass >::GetShiftValue

returns the shift value used

If type_spectrum_wanted is set to CENTERED_EIGENVALUES, we search closest eigenvalues to the shift value. Matrix (A - (shift + i shift_imag)*I)^{-1} will be used instead of A

Definition at line 274 of file EigenvalueSolver.cxx.

◆ InitMatrix() [1/2]

template<class T , class MatStiff , class MatMass >
void Seldon::EigenProblem_Base< T, MatStiff, MatMass >::InitMatrix ( MatStiff &  K)

initialization of a standard eigenvalue problem

Stiffness matrix K is given in argument. we will search (lambda, x) such as K x = lambda x

Definition at line 93 of file EigenvalueSolver.cxx.

◆ InitMatrix() [2/2]

template<class T , class MatStiff , class MatMass >
void Seldon::EigenProblem_Base< T, MatStiff, MatMass >::InitMatrix ( MatStiff &  K,
MatMass &  M 
)

initialization of a generalized eigenvalue problem

Mass matrix M and stiffness matrix K are given in argument we will search (lambda, x) such as K x = lambda M x

Definition at line 117 of file EigenvalueSolver.cxx.

◆ SetTypeSpectrum() [1/2]

template<class T , class MatStiff , class MatMass >
void Seldon::EigenProblem_Base< T, MatStiff, MatMass >::SetTypeSpectrum ( int  type,
const complex< T > &  val,
int  type_sort = SORTED_MODULUS 
)

sets which eigenvalues are searched

You can ask small eigenvalues, large, or eigenvalues close to the shift.

Definition at line 332 of file EigenvalueSolver.cxx.

◆ SetTypeSpectrum() [2/2]

template<class T , class MatStiff , class MatMass >
void Seldon::EigenProblem_Base< T, MatStiff, MatMass >::SetTypeSpectrum ( int  type,
const T &  val,
int  type_sort = SORTED_MODULUS 
)

sets which eigenvalues are searched

You can ask small eigenvalues, large, or eigenvalues close to the shift.

Definition at line 317 of file EigenvalueSolver.cxx.

Friends And Related Function Documentation

◆ ApplyScalingEigenvec [1/4]

template<class T , class MatStiff = Matrix<T>, class MatMass = Matrix<double>>
template<class EigenPb , class Vector1 , class Matrix1 , class T0 >
void ApplyScalingEigenvec ( EigenPb &  var,
Vector1 &  eigen_values,
Vector1 &  lambda_imag,
Matrix1 &  eigen_vectors,
const complex< T0 > &  shiftr,
const complex< T0 > &  shifti 
)
friend

modification of eigenvectors to take into account scaling by mass matrix

One may desire to use matrix D^-1/2 K D^-1/2 or L^-1 K L^-T instead of K in order to solve a standard eigenvalue problem instead of a generalized one. => eigenvectors are recovered by multiplying them by matrix D^1/2 or by L^T with this function

Definition at line 983 of file EigenvalueSolver.cxx.

◆ ApplyScalingEigenvec [2/4]

template<class T , class MatStiff = Matrix<T>, class MatMass = Matrix<double>>
template<class EigenPb , class Vector1 , class Matrix1 , class T0 >
void ApplyScalingEigenvec ( EigenPb &  var,
Vector1 &  eigen_values,
Vector1 &  lambda_imag,
Matrix1 &  eigen_vectors,
const T0 &  shiftr,
const T0 &  shifti 
)
friend

modification of eigenvectors to take into account scaling by mass matrix

One may desire to use matrix D^-1/2 K D^-1/2 or L^-1 K L^-T instead of K in order to solve a standard eigenvalue problem instead of a generalized one. => eigenvectors are recovered by multiplying them by matrix D^1/2 or by L^T with this function

Definition at line 846 of file EigenvalueSolver.cxx.

◆ ApplyScalingEigenvec [3/4]

template<class T , class MatStiff = Matrix<T>, class MatMass = Matrix<double>>
template<class T0 , class Prop , class Storage >
void ApplyScalingEigenvec ( EigenProblem_Base< complex< T0 > > &  var,
Vector< complex< T0 > > &  eigen_values,
Vector< complex< T0 > > &  lambda_imag,
Matrix< complex< T0 >, Prop, Storage > &  eigen_vectors,
const complex< T0 > &  shiftr,
const complex< T0 > &  shifti 
)
friend

modification of eigenvectors to take into account scaling by mass matrix

One may desire to use matrix D^-1/2 K D^-1/2 or L^-1 K L^-T instead of K in order to solve a standard eigenvalue problem instead of a generalized one. => eigenvectors are recovered by multiplying them by matrix D^1/2 or by L^T with this function

Definition at line 1618 of file VirtualEigenvalueSolver.cxx.

◆ ApplyScalingEigenvec [4/4]

template<class T , class MatStiff = Matrix<T>, class MatMass = Matrix<double>>
template<class T0 , class Prop , class Storage >
void ApplyScalingEigenvec ( EigenProblem_Base< T0 > &  var,
Vector< T0 > &  eigen_values,
Vector< T0 > &  lambda_imag,
Matrix< T0, Prop, Storage > &  eigen_vectors,
const T0 &  shiftr,
const T0 &  shifti 
)
friend

modification of eigenvectors to take into account scaling by mass matrix

One may desire to use matrix D^-1/2 K D^-1/2 or L^-1 K L^-T instead of K in order to solve a standard eigenvalue problem instead of a generalized one. => eigenvectors are recovered by multiplying them by matrix D^1/2 or by L^T with this function

Definition at line 1468 of file VirtualEigenvalueSolver.cxx.

Member Data Documentation

◆ diagonal_mass

template<class T , class MatStiff = Matrix<T>, class MatMass = Matrix<double>>
bool Seldon::EigenProblem_Base< T, MatStiff, MatMass >::diagonal_mass
protected

if true, the generalized problem is reduced to a standard problem

if M is diagonal, one can seek eigenvalues of the standard problem M^-1/2 K M^-1/2 x = lambda x

Definition at line 98 of file EigenvalueSolver.hxx.

◆ nb_add_eigenvalues

template<class T , class MatStiff = Matrix<T>, class MatMass = Matrix<double>>
int Seldon::EigenProblem_Base< T, MatStiff, MatMass >::nb_add_eigenvalues
protected

additional number of eigenvalues

Sometimes Arpack finds more converged eigenvalues than asked it is needed to store these eigenvalues and eigenvalues to avoid segmentation fault

Definition at line 77 of file EigenvalueSolver.hxx.

◆ use_cholesky

template<class T , class MatStiff = Matrix<T>, class MatMass = Matrix<double>>
bool Seldon::EigenProblem_Base< T, MatStiff, MatMass >::use_cholesky
protected

if true, the generalized problem is reduced to a standard problem

If matrix M is symmetric definite positive, one may compute Cholesky factorisation of M = L L^T, and find eigenvalues of the standard problem : L^-1 K L^-T x = lambda x

Definition at line 91 of file EigenvalueSolver.hxx.


The documentation for this class was generated from the following files: