Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
List of all members
Anasazi::LOBPCGSolMgr< ScalarType, MV, OP > Class Template Reference

User interface for the LOBPCG eigensolver. More...

#include <AnasaziLOBPCGSolMgr.hpp>

Inheritance diagram for Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >:
Anasazi::SolverManager< ScalarType, MV, OP >

Public Member Functions

Constructors/Destructor
 LOBPCGSolMgr (const Teuchos::RCP< Eigenproblem< ScalarType, MV, OP > > &problem, Teuchos::ParameterList &pl)
 Basic constructor for LOBPCGSolMgr. More...
 
virtual ~LOBPCGSolMgr ()
 Destructor. More...
 
Accessor methods
const Eigenproblem< ScalarType,
MV, OP > & 
getProblem () const
 Return the eigenvalue problem. More...
 
int getNumIters () const
 Get the iteration count for the most recent call to solve(). More...
 
Teuchos::Array< Teuchos::RCP
< Teuchos::Time > > 
getTimers () const
 Return the timers for this object. More...
 
Solver application methods
ReturnType solve ()
 This method performs possibly repeated calls to the underlying eigensolver's iterate() routine until the problem has been solved (as decided by the solver manager) or the solver manager decides to quit. More...
 
void setGlobalStatusTest (const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &global)
 Set the status test defining global convergence. More...
 
const Teuchos::RCP< StatusTest
< ScalarType, MV, OP > > & 
getGlobalStatusTest () const
 Get the status test defining global convergence. More...
 
void setLockingStatusTest (const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &locking)
 Set the status test defining locking. More...
 
const Teuchos::RCP< StatusTest
< ScalarType, MV, OP > > & 
getLockingStatusTest () const
 Get the status test defining locking. More...
 
void setDebugStatusTest (const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &debug)
 Set the status test for debugging. More...
 
const Teuchos::RCP< StatusTest
< ScalarType, MV, OP > > & 
getDebugStatusTest () const
 Get the status test for debugging. More...
 
- Public Member Functions inherited from Anasazi::SolverManager< ScalarType, MV, OP >
 SolverManager ()
 Empty constructor. More...
 
virtual ~SolverManager ()
 Destructor. More...
 

Detailed Description

template<class ScalarType, class MV, class OP>
class Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >

User interface for the LOBPCG eigensolver.

This class provides a user interface for the LOBPCG (Locally Optimal Block Preconditioned Conjugate Gradient) eigensolver. It provides the following features:

Examples:
LOBPCGCustomStatusTest.cpp, LOBPCGEpetra.cpp, LOBPCGEpetraEx.cpp, LOBPCGEpetraExGen.cpp, LOBPCGEpetraExGenPrecIfpack.cpp, and LOBPCGEpetraExGenShifted.cpp.

Definition at line 174 of file AnasaziLOBPCGSolMgr.hpp.

Constructor & Destructor Documentation

template<class ScalarType , class MV , class OP >
Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >::LOBPCGSolMgr ( const Teuchos::RCP< Eigenproblem< ScalarType, MV, OP > > &  problem,
Teuchos::ParameterList pl 
)

Basic constructor for LOBPCGSolMgr.

This constructor accepts the Eigenproblem to be solved in addition to a parameter list of options for the solver manager. These options include the following:

  • Solver parameters
    • "Which" - a string specifying the desired eigenvalues: SM, LM, SR or LR. Default: "SR"
    • "Block Size" - an int specifying the block size to be used by the underlying LOBPCG solver. Default: problem->getNEV()
    • "Full Ortho" - a bool specifying whether the underlying solver should employ the full orthogonalization scheme. Default: true
    • "Recover" - a bool specifying whether the solver manager should attempt to recover in the case of a LOBPCGRitzFailure when full orthogonalization is disabled. Default: true
    • "Verbosity" - a sum of MsgType specifying the verbosity. Default: Errors
    • "Init" - a LOBPCGState<ScalarType,MV> struct used to initialize the LOBPCG eigensolver.
    • "Output Stream" - a reference-counted pointer to the formatted output stream where all solver output is sent. Default: Teuchos::getFancyOStream ( Teuchos::rcpFromRef (std::cout) )
    • "Output Processor" - an int specifying the MPI processor that will print solver/timer details. Default: 0
  • Convergence parameters (if using default convergence test; see setGlobalStatusTest())
    • "Maximum Iterations" - an int specifying the maximum number of iterations the underlying solver is allowed to perform. Default: 100
    • "Convergence Tolerance" - a MagnitudeType specifying the level that residual norms must reach to decide convergence. Default: machine precision.
    • "Relative Convergence Tolerance" - a bool specifying whether residuals norms should be scaled by their eigenvalues for the purposing of deciding convergence. Default: true
    • "Convergence Norm" - a string specifying the norm for convergence testing: "2" or "M"
  • Locking parameters (if using default locking test; see setLockingStatusTest())
    • "Use Locking" - a bool specifying whether the algorithm should employ locking of converged eigenpairs. Default: false
    • "Max Locked" - an int specifying the maximum number of eigenpairs to be locked. Default: problem->getNEV()
    • "Locking Quorum" - an int specifying the number of eigenpairs that must meet the locking criteria before locking actually occurs. Default: 1
    • "Locking Tolerance" - a MagnitudeType specifying the level that residual norms must reach to decide locking. Default: 0.1*convergence tolerance
    • "Relative Locking Tolerance" - a bool specifying whether residuals norms should be scaled by their eigenvalues for the purposing of deciding locking. Default: true
    • "Locking Norm" - a string specifying the norm for locking testing: "2" or "M"

Definition at line 330 of file AnasaziLOBPCGSolMgr.hpp.

template<class ScalarType, class MV, class OP>
virtual Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >::~LOBPCGSolMgr ( )
inlinevirtual

Destructor.

Definition at line 219 of file AnasaziLOBPCGSolMgr.hpp.

Member Function Documentation

template<class ScalarType, class MV, class OP>
const Eigenproblem<ScalarType,MV,OP>& Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >::getProblem ( ) const
inlinevirtual

Return the eigenvalue problem.

Implements Anasazi::SolverManager< ScalarType, MV, OP >.

Definition at line 226 of file AnasaziLOBPCGSolMgr.hpp.

template<class ScalarType, class MV, class OP>
int Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >::getNumIters ( ) const
inlinevirtual

Get the iteration count for the most recent call to solve().

Implements Anasazi::SolverManager< ScalarType, MV, OP >.

Definition at line 231 of file AnasaziLOBPCGSolMgr.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::Array<Teuchos::RCP<Teuchos::Time> > Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >::getTimers ( ) const
inlinevirtual

Return the timers for this object.

The timers are ordered as follows:

  • time spent in solve() routine
  • time spent locking converged eigenvectors

Reimplemented from Anasazi::SolverManager< ScalarType, MV, OP >.

Definition at line 241 of file AnasaziLOBPCGSolMgr.hpp.

template<class ScalarType , class MV , class OP >
ReturnType Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >::solve ( )
virtual

This method performs possibly repeated calls to the underlying eigensolver's iterate() routine until the problem has been solved (as decided by the solver manager) or the solver manager decides to quit.

This method calls LOBPCG::iterate(), which will return either because a specially constructed status test evaluates to Passed or an exception is thrown.

A return from LOBPCG::iterate() signifies one of the following scenarios:

  • the maximum number of iterations has been exceeded. In this scenario, the solver manager will place
    all converged eigenpairs into the eigenproblem and return Unconverged.
  • the locking conditions have been met. In this scenario, some of the current eigenpairs will be removed
    from the eigensolver and placed into auxiliary storage. The eigensolver will be restarted with the remaining
    eigenpairs and some random information to replace the removed eigenpairs.
  • global convergence has been met. In this case, the most significant NEV eigenpairs in the solver and locked storage
    have met the convergence criterion. (Here, NEV refers to the number of eigenpairs requested by the Eigenproblem.)
    In this scenario, the solver manager will return Converged.
  • an LOBPCGRitzFailure exception has been thrown. If full orthogonalization is enabled and recovery from this exception
    is requested, the solver manager will attempt to recover from this exception by gathering the current eigenvectors,
    preconditioned residual, and search directions from the eigensolver, orthogonormalizing the basis composed of these
    three, projecting the eigenproblem, and restarting the eigensolver with the solution of the project eigenproblem. Any
    additional failure that occurs during this recovery effort will result in the eigensolver returning Unconverged.
Returns
ReturnType specifying:
  • Converged: the eigenproblem was solved to the specification required by the solver manager.
  • Unconverged: the eigenproblem was not solved to the specification desired by the solver manager

Implements Anasazi::SolverManager< ScalarType, MV, OP >.

Definition at line 475 of file AnasaziLOBPCGSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >::setGlobalStatusTest ( const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &  global)

Set the status test defining global convergence.

Definition at line 1063 of file AnasaziLOBPCGSolMgr.hpp.

template<class ScalarType , class MV , class OP >
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > & Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >::getGlobalStatusTest ( ) const

Get the status test defining global convergence.

Definition at line 1071 of file AnasaziLOBPCGSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >::setLockingStatusTest ( const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &  locking)

Set the status test defining locking.

Definition at line 1093 of file AnasaziLOBPCGSolMgr.hpp.

template<class ScalarType , class MV , class OP >
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > & Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >::getLockingStatusTest ( ) const

Get the status test defining locking.

Definition at line 1101 of file AnasaziLOBPCGSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >::setDebugStatusTest ( const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &  debug)

Set the status test for debugging.

Definition at line 1078 of file AnasaziLOBPCGSolMgr.hpp.

template<class ScalarType , class MV , class OP >
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > & Anasazi::LOBPCGSolMgr< ScalarType, MV, OP >::getDebugStatusTest ( ) const

Get the status test for debugging.

Definition at line 1086 of file AnasaziLOBPCGSolMgr.hpp.


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