Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
Belos::BlockGCRODRSolMgr< ScalarType, MV, OP > Class Template Reference

A solver manager for the Block GCRO-DR (Block Recycling GMRES) linear solver. More...

#include <BelosBlockGCRODRSolMgr.hpp>

Inheritance diagram for Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >:
Inheritance graph
[legend]

Private Types

typedef MultiVecTraits
< ScalarType, MV > 
MVT
 
typedef OperatorTraits
< ScalarType, MV, OP > 
OPT
 
typedef Teuchos::ScalarTraits
< ScalarType > 
SCT
 
typedef Teuchos::ScalarTraits
< ScalarType >::magnitudeType 
MagnitudeType
 
typedef Teuchos::ScalarTraits
< MagnitudeType
MT
 
typedef Teuchos::ScalarTraits
< MagnitudeType
SMT
 
typedef OrthoManagerFactory
< ScalarType, MV, OP > 
ortho_factory_type
 
typedef
Teuchos::SerialDenseMatrix
< int, ScalarType > 
SDM
 
typedef
Teuchos::SerialDenseVector
< int, ScalarType > 
SDV
 

Private Member Functions

void init ()
 
void initializeStateStorage ()
 
void buildRecycleSpaceKryl (int &keff, Teuchos::RCP< BlockGmresIter< ScalarType, MV, OP > > block_gmres_iter)
 
void buildRecycleSpaceAugKryl (Teuchos::RCP< BlockGCRODRIter< ScalarType, MV, OP > > gcrodr_iter)
 
int getHarmonicVecsKryl (int m, const SDM &HH, SDM &PP)
 
int getHarmonicVecsAugKryl (int keff, int m, const SDM &HH, const Teuchos::RCP< const MV > &VV, SDM &PP)
 
void sort (std::vector< MagnitudeType > &dlist, int n, std::vector< int > &iperm)
 

Private Attributes

Teuchos::LAPACK< int, ScalarType > lapack
 
Teuchos::RCP< LinearProblem
< ScalarType, MV, OP > > 
problem_
 The current linear problem to solve. More...
 
Teuchos::RCP< OutputManager
< ScalarType > > 
printer_
 
Teuchos::RCP< std::ostream > outputStream_
 
Teuchos::RCP< StatusTest
< ScalarType, MV, OP > > 
sTest_
 
Teuchos::RCP
< StatusTestMaxIters
< ScalarType, MV, OP > > 
maxIterTest_
 
Teuchos::RCP< StatusTest
< ScalarType, MV, OP > > 
convTest_
 
Teuchos::RCP
< StatusTestGenResNorm
< ScalarType, MV, OP > > 
expConvTest_
 
Teuchos::RCP
< StatusTestGenResNorm
< ScalarType, MV, OP > > 
impConvTest_
 
Teuchos::RCP< StatusTestOutput
< ScalarType, MV, OP > > 
outputTest_
 
ortho_factory_type orthoFactory_
 Factory for creating MatOrthoManager subclass instances. More...
 
Teuchos::RCP< MatOrthoManager
< ScalarType, MV, OP > > 
ortho_
 Orthogonalization manager. More...
 
Teuchos::RCP
< Teuchos::ParameterList
params_
 This solver's current parameter list. More...
 
Teuchos::RCP< const
Teuchos::ParameterList
defaultParams_
 Default parameter list. More...
 
MagnitudeType convTol_
 
MagnitudeType orthoKappa_
 
MagnitudeType achievedTol_
 
int blockSize_
 
int maxRestarts_
 
int maxIters_
 
int numIters_
 
int verbosity_
 
int outputStyle_
 
int outputFreq_
 
bool adaptiveBlockSize_
 
std::string orthoType_
 
std::string recycleMethod_
 
std::string impResScale_
 
std::string expResScale_
 
std::string label_
 
int numBlocks_
 
int recycledBlocks_
 
int keff
 
Teuchos::RCP< MV > R_
 
Teuchos::RCP< MV > V_
 
Teuchos::RCP< MV > U_
 
Teuchos::RCP< MV > C_
 
Teuchos::RCP< MV > U1_
 
Teuchos::RCP< MV > C1_
 
Teuchos::RCP< SDMG_
 
Teuchos::RCP< SDMH_
 
Teuchos::RCP< SDMB_
 
Teuchos::RCP< SDMPP_
 
Teuchos::RCP< SDMHP_
 
std::vector< ScalarType > tau_
 
std::vector< ScalarType > work_
 
Teuchos::RCP< SDMF_
 
std::vector< int > ipiv_
 
Teuchos::RCP< Teuchos::TimetimerSolve_
 Timer for solve(). More...
 
bool isSet_
 Whether setParameters() successfully finished setting parameters. More...
 
bool loaDetected_
 Whether a loss of accuracy was detected during the solve. More...
 
bool builtRecycleSpace_
 Whether we have generated or regenerated a recycle space yet this solve. More...
 

Static Private Attributes

static const bool adaptiveBlockSize_default_ = true
 
static const std::string recycleMethod_default_ = "harmvecs"
 

Constructors/Destructor

 BlockGCRODRSolMgr ()
 Default constructor. More...
 
 BlockGCRODRSolMgr (const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl)
 Basic constructor for GCRODRSolMgr. More...
 
virtual ~BlockGCRODRSolMgr ()
 Destructor. More...
 

Implementation of the Teuchos::Describable interface

std::string description () const
 A description of the Block GCRODR solver manager. More...
 

Accessor methods

const LinearProblem
< ScalarType, MV, OP > & 
getProblem () const
 Get current linear problem being solved for in this object. More...
 
Teuchos::RCP< const
Teuchos::ParameterList
getValidParameters () const
 Get a parameter list containing the valid parameters for this object. More...
 
Teuchos::RCP< const
Teuchos::ParameterList
getCurrentParameters () const
 Get a parameter list containing the current parameters for this object. More...
 
MagnitudeType achievedTol () const
 Get the residual for the most recent call to solve(). More...
 
int getNumIters () const
 Get the iteration count for the most recent call to solve(). More...
 
bool isLOADetected () const
 Whether a loss of accuracy was detected during the most recent solve. More...
 

Set methods

void setProblem (const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)
 Set the linear problem to solve on the next call to solve(). More...
 
void setParameters (const Teuchos::RCP< Teuchos::ParameterList > &params)
 Set the parameters the solver should use to solve the linear problem. More...
 

Reset methods

void reset (const ResetType type)
 Performs a reset of the solver manager specified by the ResetType. More...
 

Solver application methods

ReturnType solve ()
 Solve the current linear problem. More...
 

Additional Inherited Members

- Public Member Functions inherited from Belos::SolverManager< ScalarType, MV, OP >
 SolverManager ()
 Empty constructor. More...
 
virtual ~SolverManager ()
 Destructor. More...
 
virtual Teuchos::RCP
< SolverManager< ScalarType,
MV, OP > > 
clone () const =0
 clone the solver manager. More...
 
virtual void setUserConvStatusTest (const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &, const typename StatusTestCombo< ScalarType, MV, OP >::ComboType &=StatusTestCombo< ScalarType, MV, OP >::SEQ)
 Set user-defined convergence status test. More...
 
virtual void setDebugStatusTest (const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &)
 Set user-defined debug status test. More...
 
- Public Member Functions inherited from Teuchos::Describable
DescribableStreamManipulatorState describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default)
 
std::ostream & operator<< (std::ostream &os, const DescribableStreamManipulatorState &d)
 
virtual void describe (FancyOStream &out, const EVerbosityLevel verbLevel=verbLevel_default) const
 
void describe (std::ostream &out, const EVerbosityLevel verbLevel=verbLevel_default) const
 
virtual ~Describable ()
 
 LabeledObject ()
 
virtual ~LabeledObject ()
 
virtual void setObjectLabel (const std::string &objectLabel)
 
virtual std::string getObjectLabel () const
 
- Static Public Attributes inherited from Teuchos::Describable
static const EVerbosityLevel verbLevel_default
 

Detailed Description

template<class ScalarType, class MV, class OP>
class Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >

A solver manager for the Block GCRO-DR (Block Recycling GMRES) linear solver.

Author
Kirk M. Soodhalter and Michael L. Parks

GCRO-DR (also called Recycling GMRES) is a variant of GMRES that can more efficiently solve sequences linear systems. It does so by "recycling" Krylov basis information from previous solves.

The original GCRO-DR algorithm can only solve one right-hand side at a time. Block GCRO-DR extends GCRO-DR so that it can solve multiple right-hand sides at a time; thus, it can solve sequences of block systems.

Definition at line 127 of file BelosBlockGCRODRSolMgr.hpp.

Member Typedef Documentation

template<class ScalarType , class MV , class OP >
typedef MultiVecTraits<ScalarType,MV> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::MVT
private

Definition at line 130 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
typedef OperatorTraits<ScalarType,MV,OP> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::OPT
private

Definition at line 131 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
typedef Teuchos::ScalarTraits<ScalarType> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::SCT
private

Definition at line 132 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
typedef Teuchos::ScalarTraits<ScalarType>::magnitudeType Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::MagnitudeType
private

Definition at line 133 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
typedef Teuchos::ScalarTraits<MagnitudeType> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::MT
private

Definition at line 134 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
typedef Teuchos::ScalarTraits<MagnitudeType> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::SMT
private

Definition at line 135 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
typedef OrthoManagerFactory<ScalarType, MV, OP> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::ortho_factory_type
private

Definition at line 136 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
typedef Teuchos::SerialDenseMatrix<int,ScalarType> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::SDM
private

Definition at line 137 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
typedef Teuchos::SerialDenseVector<int,ScalarType> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::SDV
private

Definition at line 138 of file BelosBlockGCRODRSolMgr.hpp.

Constructor & Destructor Documentation

template<class ScalarType , class MV , class OP >
Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::BlockGCRODRSolMgr ( )

Default constructor.

This constructor sets up the solver with default parameters. The linear problem must be passed in using setProblem() before solve() is called on this object. The solver values can be changed using setParameters().

Definition at line 438 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::BlockGCRODRSolMgr ( const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &  problem,
const Teuchos::RCP< Teuchos::ParameterList > &  pl 
)

Basic constructor for GCRODRSolMgr.

This constructor accepts the LinearProblem to be solved in addition to a parameter list of options for the solver manager. Some of the more important options include the following:

  • "Num Blocks": an int specifying the number of blocks allocated for the Krylov basis. Default: 50.
  • "Block Size": an int specifying the number of right hand sides being solved at a time.
  • "Num Recycled Blocks": an int specifying the number of blocks allocated for the Krylov basis. Default: 5.
  • "Maximum Iterations": an int specifying the maximum number of iterations the underlying solver is allowed to perform. Default: 5000.
  • "Maximum Restarts": an int specifying the maximum number of restarts the underlying solver is allowed to perform. Default: 100.
  • "Orthogonalization": an std::string specifying the desired orthogonalization. Currently supported values: "DGKS", "ICGS", "IMGS", and "TSQR" (if Belos was built with TSQR support). Default: "ICGS".
  • "Orthogonalization Parameters": a ParameterList or RCP<(const) ParameterList> of parameters specific to the type of orthogonalization used. Defaults are set automatically.
  • "Verbosity": a sum of MsgType specifying the verbosity. Default: Belos::Errors.
  • "Output Style": a OutputType specifying the style of output. Default: Belos::General.
  • "Convergence Tolerance": a MagnitudeType specifying the level that residual norms must reach to decide convergence. Default: 1e-8.

Other supported options:

  • "Output Frequency": an int specifying how often (in terms of number of iterations) convergence information should be output to the output stream. Default: -1 (never output convergence information).
  • "Output Stream": a reference-counted pointer to the output stream where all solver output is sent. Default stream is std::cout (stdout, in C terms). For stderr, supply Teuchos::rcp(&std::cerr, false).
  • "Implicit Residual Scaling": the type of scaling used in the implicit residual convergence test. Default: "Norm of Preconditioned Initial Residual".
  • "Explicit Residual Scaling": the type of scaling used in the explicit residual convergence test. Default: "Norm of Initial Residual".
  • "Timer Label": the string to use as a prefix for the timer labels. Default: "Belos"
  • "Orthogonalization Constant": a MagnitudeType corresponding to the "depTol" parameter of DGKS orthogonalization. Ignored unless DGKS orthogonalization is used. DGKS decides the default value.

Definition at line 445 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
virtual Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::~BlockGCRODRSolMgr ( )
inlinevirtual

Destructor.

Definition at line 180 of file BelosBlockGCRODRSolMgr.hpp.

Member Function Documentation

template<class ScalarType , class MV , class OP >
std::string Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::description ( ) const
virtual

A description of the Block GCRODR solver manager.

Reimplemented from Teuchos::Describable.

Definition at line 541 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
const LinearProblem<ScalarType,MV,OP>& Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::getProblem ( ) const
inlinevirtual

Get current linear problem being solved for in this object.

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

Definition at line 196 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP< const Teuchos::ParameterList > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::getValidParameters ( ) const
virtual

Get a parameter list containing the valid parameters for this object.

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

Definition at line 556 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<const Teuchos::ParameterList> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::getCurrentParameters ( ) const
inlinevirtual

Get a parameter list containing the current parameters for this object.

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

Definition at line 204 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::achievedTol ( ) const
inlinevirtual

Get the residual for the most recent call to solve().

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

Definition at line 207 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::getNumIters ( ) const
inlinevirtual

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

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

Definition at line 212 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::isLOADetected ( ) const
inlinevirtual

Whether a loss of accuracy was detected during the most recent solve.

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

Definition at line 215 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::setProblem ( const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &  problem)
inlinevirtual

Set the linear problem to solve on the next call to solve().

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

Definition at line 223 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::setParameters ( const Teuchos::RCP< Teuchos::ParameterList > &  params)
virtual

Set the parameters the solver should use to solve the linear problem.

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

Definition at line 643 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::reset ( const ResetType  type)
inlinevirtual

Performs a reset of the solver manager specified by the ResetType.

This informs the solver manager that the solver should prepare for the next call to solve by resetting certain elements of the iterative solver strategy.

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

Definition at line 252 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
ReturnType Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::solve ( )
virtual

Solve the current linear problem.

This method performs possibly repeated calls to the underlying linear solver'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 BlockGCRODRIter::iterate(), which will return either because a specially constructed status test evaluates to Passed or an exception is thrown.

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

  • the maximum number of restarts has been exceeded. In this scenario, the current solutions to the linear system will be placed in the linear problem and return Unconverged.
  • global convergence has been met. In this case, the current solutions to the linear system will be placed in the linear problem and the solver manager will return Converged.- Converged: the linear problem was solved to the specification required by the solver manager.
  • Unconverged: the linear problem was not solved to the specification desired by the solver manager.

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

Definition at line 1801 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::init ( )
private

Definition at line 464 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::initializeStateStorage ( )
private

Definition at line 1058 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::buildRecycleSpaceKryl ( int &  keff,
Teuchos::RCP< BlockGmresIter< ScalarType, MV, OP > >  block_gmres_iter 
)
private

Definition at line 1210 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::buildRecycleSpaceAugKryl ( Teuchos::RCP< BlockGCRODRIter< ScalarType, MV, OP > >  gcrodr_iter)
private

Definition at line 1304 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::getHarmonicVecsKryl ( int  m,
const SDM HH,
SDM PP 
)
private

Definition at line 1606 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::getHarmonicVecsAugKryl ( int  keff,
int  m,
const SDM HH,
const Teuchos::RCP< const MV > &  VV,
SDM PP 
)
private

Definition at line 1473 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::sort ( std::vector< MagnitudeType > &  dlist,
int  n,
std::vector< int > &  iperm 
)
private

Definition at line 1741 of file BelosBlockGCRODRSolMgr.hpp.

Member Data Documentation

template<class ScalarType , class MV , class OP >
Teuchos::LAPACK<int,ScalarType> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::lapack
private

Definition at line 325 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::problem_
private

The current linear problem to solve.

Definition at line 328 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<OutputManager<ScalarType> > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::printer_
private

Definition at line 331 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<std::ostream> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::outputStream_
private

Definition at line 332 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTest<ScalarType,MV,OP> > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::sTest_
private

Definition at line 335 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::maxIterTest_
private

Definition at line 336 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTest<ScalarType,MV,OP> > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::convTest_
private

Definition at line 337 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::expConvTest_
private

Definition at line 338 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::impConvTest_
private

Definition at line 338 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::outputTest_
private

Definition at line 339 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
ortho_factory_type Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::orthoFactory_
private

Factory for creating MatOrthoManager subclass instances.

Definition at line 342 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::ortho_
private

Orthogonalization manager.

This is created by the OrthoManagerFactory instance. The pointer may be invalidated if this solver's parameters are changed.

Definition at line 349 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<Teuchos::ParameterList> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::params_
private

This solver's current parameter list.

Definition at line 352 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<const Teuchos::ParameterList> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::defaultParams_
mutableprivate

Default parameter list.

This is declared "mutable" so that it can be created on demand by getValidParameters(). The SolverManager interface forces getValidParameters() to be a const method.

Using the "mutable" keyword is ugly, but the previous caching solution involved static method data in getValidParameters(). That solution made creating multiple instances of this solver in different threads not safe.

Definition at line 364 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
const bool Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::adaptiveBlockSize_default_ = true
staticprivate

Definition at line 367 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
const std::string Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::recycleMethod_default_ = "harmvecs"
staticprivate

Definition at line 368 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::convTol_
private

Definition at line 371 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::orthoKappa_
private

Definition at line 371 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::achievedTol_
private

Definition at line 371 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::blockSize_
private

Definition at line 372 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::maxRestarts_
private

Definition at line 372 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::maxIters_
private

Definition at line 372 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::numIters_
private

Definition at line 372 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::verbosity_
private

Definition at line 373 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::outputStyle_
private

Definition at line 373 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::outputFreq_
private

Definition at line 373 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::adaptiveBlockSize_
private

Definition at line 374 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
std::string Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::orthoType_
private

Definition at line 375 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
std::string Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::recycleMethod_
private

Definition at line 375 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
std::string Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::impResScale_
private

Definition at line 376 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
std::string Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::expResScale_
private

Definition at line 376 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
std::string Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::label_
private

Definition at line 377 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::numBlocks_
private

Definition at line 384 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::recycledBlocks_
private

Definition at line 384 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::keff
private

Definition at line 386 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<MV> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::R_
private

Definition at line 389 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<MV> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::V_
private

Definition at line 392 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<MV> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::U_
private

Definition at line 395 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<MV> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::C_
private

Definition at line 395 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<MV> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::U1_
private

Definition at line 398 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<MV> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::C1_
private

Definition at line 398 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<SDM > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::G_
private

Definition at line 401 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<SDM > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::H_
private

Definition at line 402 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<SDM > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::B_
private

Definition at line 403 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<SDM > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::PP_
private

Definition at line 404 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<SDM > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::HP_
private

Definition at line 405 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
std::vector<ScalarType> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::tau_
private

Definition at line 406 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
std::vector<ScalarType> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::work_
private

Definition at line 407 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<SDM > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::F_
private

Definition at line 408 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
std::vector<int> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::ipiv_
private

Definition at line 409 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<Teuchos::Time> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::timerSolve_
private

Timer for solve().

Definition at line 412 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::isSet_
private

Whether setParameters() successfully finished setting parameters.

Definition at line 415 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::loaDetected_
private

Whether a loss of accuracy was detected during the solve.

Definition at line 418 of file BelosBlockGCRODRSolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::builtRecycleSpace_
private

Whether we have generated or regenerated a recycle space yet this solve.

Definition at line 421 of file BelosBlockGCRODRSolMgr.hpp.


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