#include <jdbsym.hh>
Public Member Functions | |
JdbSym (concepts::Operator< Real > &stiff, concepts::Operator< Real > &mass, Real tol, uint maxit=150, Real tau=0.0, uint jdtype=1, uint kmax=1, concepts::SolverFabric< Real > *fabric=0, const concepts::Array< concepts::Vector< Real > * > *start=0) | |
virtual const concepts::Array< Real > & | getEV () |
virtual const concepts::Array< concepts::Vector< Real > * > & | getEF () |
virtual uint | iterations () const |
Returns the number of iterations. | |
virtual uint | converged () const |
Returns the number of converged eigen pairs. | |
Protected Member Functions | |
virtual std::ostream & | info (std::ostream &os) const |
Returns information in an output stream. | |
Eigenvalue solver using JDBSYM.
JDBSYM is an implementation of the Jacobi-Davidson method optimized for symmetric eigenvalue problems. It solves eigenproblems of the form and with or without preconditioning, where A is symmetric and B is symmetric positive definite. The implementation supports blocking.
eigensolver::JdbSym::JdbSym | ( | concepts::Operator< Real > & | stiff, |
concepts::Operator< Real > & | mass, | ||
Real | tol, | ||
uint | maxit = 150 , |
||
Real | tau = 0.0 , |
||
uint | jdtype = 1 , |
||
uint | kmax = 1 , |
||
concepts::SolverFabric< Real > * | fabric = 0 , |
||
const concepts::Array< concepts::Vector< Real > * > * | start = 0 |
||
) |
Constructor. Initializes the JdbSym class. As parameter it takes the few most important variables that are needed for the Jacobi-Davidson algorithm. The rest of the parameters are given as default and work in most cases. The class then solves the problem where A, B are the stiffness and the mass matrices respectively.
The computation is only perfomed if getEV
or getEF
is called.
stiff | Stiffness matrix A |
mass | Mass matrix B |
fabric | Solver fabric for a linear solver (preconditioner for the shifted stiffness matrix) |
tol | Convergence tolerance for the eigenpairs. For a pair convergence is defined by |
maxit | Maximal number of iterations |
tau | Target value of Jacobi-Davidson algorithm. The code will find the kmax eigenvalues closest to tau . |
jdtype | Type of solver required. An older solver and a newer solver are possible to use. |
kmax | Number of eigenpairs to be computed |
start | Starting vectors. Used to build the initial search subspace |
|
inlinevirtual |
Returns the number of converged eigen pairs.
Implements eigensolver::EigenSolver< Real >.
|
virtual |
Implements eigensolver::EigenSolver< Real >.
|
virtual |
Returns an array with the eigen values
Implements eigensolver::EigenSolver< Real >.
|
protectedvirtual |
Returns information in an output stream.
Reimplemented from eigensolver::EigenSolver< Real >.
|
inlinevirtual |
Returns the number of iterations.
Implements eigensolver::EigenSolver< Real >.