45 #ifndef THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP 
   46 #define THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP 
   50 #include "Thyra_BelosLinearOpWithSolve.hpp" 
   51 #include "Thyra_ScaledAdjointLinearOpBase.hpp" 
   53 #include "BelosBlockGmresSolMgr.hpp" 
   54 #include "BelosPseudoBlockGmresSolMgr.hpp" 
   55 #include "BelosBlockCGSolMgr.hpp" 
   56 #include "BelosPseudoBlockCGSolMgr.hpp" 
   57 #include "BelosPseudoBlockStochasticCGSolMgr.hpp" 
   58 #include "BelosGCRODRSolMgr.hpp" 
   59 #include "BelosRCGSolMgr.hpp" 
   60 #include "BelosMinresSolMgr.hpp" 
   61 #include "BelosTFQMRSolMgr.hpp" 
   77 template<
class Scalar>
 
   78 const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverType_name = 
"Solver Type";
 
   79 template<
class Scalar>
 
   80 const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverType_default = 
"Pseudo Block GMRES";
 
   81 template<
class Scalar>
 
   82 const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverTypes_name = 
"Solver Types";
 
   83 template<
class Scalar>
 
   84 const std::string BelosLinearOpWithSolveFactory<Scalar>::BlockGMRES_name = 
"Block GMRES";
 
   85 template<
class Scalar>
 
   86 const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockGMRES_name = 
"Pseudo Block GMRES";
 
   87 template<
class Scalar>
 
   88 const std::string BelosLinearOpWithSolveFactory<Scalar>::BlockCG_name = 
"Block CG";
 
   89 template<
class Scalar>
 
   90 const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockCG_name = 
"Pseudo Block CG";
 
   91 template<
class Scalar>
 
   92 const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockStochasticCG_name = 
"Pseudo Block Stochastic CG";
 
   93 template<
class Scalar>
 
   94 const std::string BelosLinearOpWithSolveFactory<Scalar>::GCRODR_name = 
"GCRODR";
 
   95 template<
class Scalar>
 
   96 const std::string BelosLinearOpWithSolveFactory<Scalar>::RCG_name = 
"RCG";
 
   97 template<
class Scalar>
 
   98 const std::string BelosLinearOpWithSolveFactory<Scalar>::MINRES_name = 
"MINRES";
 
   99 template<
class Scalar>
 
  100 const std::string BelosLinearOpWithSolveFactory<Scalar>::TFQMR_name = 
"TFQMR";
 
  101 template<
class Scalar>
 
  102 const std::string BelosLinearOpWithSolveFactory<Scalar>::ConvergenceTestFrequency_name = 
"Convergence Test Frequency";
 
  105 const std::string LeftPreconditionerIfUnspecified_name = 
"Left Preconditioner If Unspecified";
 
  111 template<
class Scalar>
 
  113   :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES),
 
  114    convergenceTestFrequency_(1)
 
  120 template<
class Scalar>
 
  122   const RCP<PreconditionerFactoryBase<Scalar> > &precFactory
 
  124   :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES)
 
  133 template<
class Scalar>
 
  140 template<
class Scalar>
 
  142   const RCP<PreconditionerFactoryBase<Scalar> > &precFactory,
 
  143   const std::string &precFactoryName
 
  148     precFactoryValidPL = precFactory->getValidParameters();
 
  149   const std::string _precFactoryName =
 
  150     ( precFactoryName != 
"" 
  152       : ( precFactoryValidPL.get() ? precFactoryValidPL->name() : 
"GENERIC PRECONDITIONER FACTORY" )
 
  154   precFactory_ = precFactory;
 
  155   precFactoryName_ = _precFactoryName;
 
  156   updateThisValidParamList();
 
  160 template<
class Scalar>
 
  168 template<
class Scalar>
 
  170   RCP<PreconditionerFactoryBase<Scalar> > *precFactory,
 
  171   std::string *precFactoryName
 
  174   if(precFactory) *precFactory = precFactory_;
 
  175   if(precFactoryName) *precFactoryName = precFactoryName_;
 
  177   precFactoryName_ = 
"";
 
  178   updateThisValidParamList();
 
  182 template<
class Scalar>
 
  184   const LinearOpSourceBase<Scalar> &fwdOpSrc
 
  187   if(precFactory_.get())
 
  188     return precFactory_->isCompatible(fwdOpSrc);
 
  193 template<
class Scalar>
 
  201 template<
class Scalar>
 
  203   const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
 
  204   LinearOpWithSolveBase<Scalar> *Op,
 
  205   const ESupportSolveUse supportSolveUse
 
  209   initializeOpImpl(fwdOpSrc,null,null,
false,Op,supportSolveUse);
 
  213 template<
class Scalar>
 
  215   const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
 
  216   LinearOpWithSolveBase<Scalar> *Op
 
  220   initializeOpImpl(fwdOpSrc,null,null,
true,Op,SUPPORT_SOLVE_UNSPECIFIED);
 
  224 template<
class Scalar>
 
  226   const EPreconditionerInputType precOpType
 
  229   if(precFactory_.get())
 
  231   return (precOpType==PRECONDITIONER_INPUT_TYPE_AS_OPERATOR);
 
  235 template<
class Scalar>
 
  237   const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
 
  238   const RCP<
const PreconditionerBase<Scalar> > &prec,
 
  239   LinearOpWithSolveBase<Scalar> *Op,
 
  240   const ESupportSolveUse supportSolveUse
 
  244   initializeOpImpl(fwdOpSrc,null,prec,
false,Op,supportSolveUse);
 
  248 template<
class Scalar>
 
  250   const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
 
  251   const RCP<
const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
 
  252   LinearOpWithSolveBase<Scalar> *Op,
 
  253   const ESupportSolveUse supportSolveUse
 
  257   initializeOpImpl(fwdOpSrc,approxFwdOpSrc,null,
false,Op,supportSolveUse);
 
  261 template<
class Scalar>
 
  263   LinearOpWithSolveBase<Scalar> *Op,
 
  264   RCP<
const LinearOpSourceBase<Scalar> > *fwdOpSrc,
 
  265   RCP<
const PreconditionerBase<Scalar> > *prec,
 
  266   RCP<
const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
 
  267   ESupportSolveUse *supportSolveUse
 
  286   if(fwdOpSrc) *fwdOpSrc = _fwdOpSrc;
 
  287   if(prec) *prec = _prec;
 
  288   if(approxFwdOpSrc) *approxFwdOpSrc = _approxFwdOpSrc;
 
  289   if(supportSolveUse) *supportSolveUse = _supportSolveUse;
 
  296 template<
class Scalar>
 
  303   paramList_ = paramList;
 
  305     Teuchos::getIntegralValue<EBelosSolverType>(*paramList_, SolverType_name);
 
  306   convergenceTestFrequency_ =
 
  307     Teuchos::getParameter<int>(*paramList_, ConvergenceTestFrequency_name);
 
  308   Teuchos::readVerboseObjectSublist(&*paramList_,
this);
 
  312 template<
class Scalar>
 
  320 template<
class Scalar>
 
  330 template<
class Scalar>
 
  338 template<
class Scalar>
 
  342   return thisValidParamList_;
 
  349 template<
class Scalar>
 
  352   std::ostringstream oss;
 
  353   oss << 
"Thyra::BelosLinearOpWithSolveFactory";
 
  364 template<
class Scalar>
 
  369   using Teuchos::tuple;
 
  370   using Teuchos::setStringToIntegralParameter;
 
  376     EBelosSolverType> >::getDummyObject());
 
  378   typedef MultiVectorBase<Scalar> MV_t;
 
  379   typedef LinearOpBase<Scalar> LO_t;
 
  381   if(validParamList.
get()==NULL) {
 
  383     setStringToIntegralParameter<EBelosSolverType>(
 
  384       SolverType_name, SolverType_default,
 
  385       "Type of linear solver algorithm to use.",
 
  388         "Pseudo Block GMRES",
 
  391         "Pseudo Block Stochastic CG",
 
  398         "Block GMRES solver for nonsymmetric linear systems.  It can also solve " 
  399         "single right-hand side systems, and can also perform Flexible GMRES " 
  400         "(where the preconditioner may change at every iteration, for example " 
  401         "for inner-outer iterations), by setting options in the \"Block GMRES\" " 
  404         "GMRES solver for nonsymmetric linear systems, that performs single " 
  405         "right-hand side solves on multiple right-hand sides at once.  It " 
  406         "exploits operator multivector multiplication in order to amortize " 
  407         "global communication costs.  Individual linear systems are deflated " 
  408         "out as they are solved.",
 
  410         "Block CG solver for symmetric (Hermitian in complex arithmetic) " 
  411         "positive definite linear systems.  It can also solve single " 
  412         "right-hand-side systems.",
 
  414         "CG solver that performs single right-hand side CG on multiple right-hand " 
  415         "sides at once.  It exploits operator multivector multiplication in order " 
  416         "to amortize global communication costs.  Individual linear systems are " 
  417         "deflated out as they are solved.",
 
  419         "stochastic CG solver that performs single right-hand side CG on multiple right-hand " 
  420         "sides at once.  It exploits operator multivector multiplication in order " 
  421         "to amortize global communication costs.  Individual linear systems are " 
  422         "deflated out as they are solved. [EXPERIMENTAL]",
 
  424         "Variant of GMRES that performs subspace recycling to accelerate " 
  425         "convergence for sequences of solves with related linear systems.  " 
  426   "Individual linear systems are deflated out as they are solved.  " 
  427   "The current implementation only supports real-valued Scalar types.",
 
  429         "CG solver for symmetric (Hermitian in complex arithmetic) positive " 
  430         "definite linear systems, that performs subspace recycling to " 
  431         "accelerate convergence for sequences of related linear systems.",
 
  433         "MINRES solver for symmetric indefinite linear systems.  It performs " 
  434         "single-right-hand-side solves on multiple right-hand sides sequentially.",
 
  436         "TFQMR (Transpose-Free QMR) solver for nonsymmetric linear systems." 
  438       tuple<EBelosSolverType>(
 
  439         SOLVER_TYPE_BLOCK_GMRES,
 
  440         SOLVER_TYPE_PSEUDO_BLOCK_GMRES,
 
  441         SOLVER_TYPE_BLOCK_CG,
 
  442         SOLVER_TYPE_PSEUDO_BLOCK_CG,
 
  443         SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG,
 
  451     validParamList->
set(ConvergenceTestFrequency_name, as<int>(1),
 
  452       "Number of linear solver iterations to skip between applying" 
  453       " user-defined convergence test.");
 
  455       LeftPreconditionerIfUnspecified_name, 
false,
 
  456       "If the preconditioner does not specify if it is left or right, and this\n" 
  457       "option is set to true, put the preconditioner on the left side.\n" 
  458       "Historically, preconditioning is on the right. Some solvers may not\n" 
  459       "support left preconditioning.");
 
  461       &solverTypesSL = validParamList->
sublist(SolverTypes_name);
 
  463       Belos::BlockGmresSolMgr<Scalar,MV_t,LO_t> mgr;
 
  464       solverTypesSL.
sublist(BlockGMRES_name).setParameters(
 
  465         *mgr.getValidParameters()
 
  469       Belos::PseudoBlockGmresSolMgr<Scalar,MV_t,LO_t> mgr;
 
  470       solverTypesSL.
sublist(PseudoBlockGMRES_name).setParameters(
 
  471         *mgr.getValidParameters()
 
  475       Belos::BlockCGSolMgr<Scalar,MV_t,LO_t> mgr;
 
  476       solverTypesSL.
sublist(BlockCG_name).setParameters(
 
  477         *mgr.getValidParameters()
 
  481       Belos::PseudoBlockCGSolMgr<Scalar,MV_t,LO_t> mgr;
 
  482       solverTypesSL.
sublist(PseudoBlockCG_name).setParameters(
 
  483         *mgr.getValidParameters()
 
  487       Belos::PseudoBlockStochasticCGSolMgr<Scalar,MV_t,LO_t> mgr;
 
  488       solverTypesSL.
sublist(PseudoBlockStochasticCG_name).setParameters(
 
  489         *mgr.getValidParameters()
 
  493       Belos::GCRODRSolMgr<Scalar,MV_t,LO_t> mgr;
 
  494       solverTypesSL.
sublist(GCRODR_name).setParameters(
 
  495         *mgr.getValidParameters()
 
  499       Belos::RCGSolMgr<Scalar,MV_t,LO_t> mgr;
 
  500       solverTypesSL.
sublist(RCG_name).setParameters(
 
  501         *mgr.getValidParameters()
 
  505       Belos::MinresSolMgr<Scalar,MV_t,LO_t> mgr;
 
  506       solverTypesSL.
sublist(MINRES_name).setParameters(
 
  507         *mgr.getValidParameters()
 
  511       Belos::TFQMRSolMgr<Scalar,MV_t,LO_t> mgr;
 
  512       solverTypesSL.
sublist(TFQMR_name).setParameters(
 
  513         *mgr.getValidParameters()
 
  517   return validParamList;
 
  521 template<
class Scalar>
 
  527   Teuchos::setupVerboseObjectSublist(&*thisValidParamList_);
 
  531 template<
class Scalar>
 
  533   const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
 
  534   const RCP<
const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
 
  535   const RCP<
const PreconditionerBase<Scalar> > &prec_in,
 
  536   const bool reusePrec,
 
  537   LinearOpWithSolveBase<Scalar> *Op,
 
  538   const ESupportSolveUse supportSolveUse
 
  543   using Teuchos::set_extra_data;
 
  545   typedef MultiVectorBase<Scalar> MV_t;
 
  546   typedef LinearOpBase<Scalar> LO_t;
 
  552     *out << 
"\nEntering Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<
">::initializeOpImpl(...) ...\n";
 
  563     fwdOp = fwdOpSrc->getOp(),
 
  564     approxFwdOp = ( approxFwdOpSrc.get() ? approxFwdOpSrc->getOp() : 
Teuchos::null );
 
  585     if(precFactory_.get()) {
 
  588           ? Teuchos::rcp_const_cast<PreconditionerBase<Scalar> >(belosOp->
extract_prec())
 
  591       bool hasExistingPrec = 
false;
 
  593         hasExistingPrec = 
true;
 
  598         hasExistingPrec = 
false;
 
  599         myPrec = precFactory_->createPrec();
 
  601       if( hasExistingPrec && reusePrec ) {
 
  606         if(approxFwdOp.get())
 
  607           precFactory_->initializePrec(approxFwdOpSrc,&*myPrec);
 
  609           precFactory_->initializePrec(fwdOpSrc,&*myPrec);
 
  619   bool oldIsExternalPrec = 
false;
 
  624   ESupportSolveUse oldSupportSolveUse = SUPPORT_SOLVE_UNSPECIFIED;
 
  626   belosOp->
uninitialize( &oldLP, NULL, &oldIterSolver, &oldFwdOpSrc,
 
  627     NULL, &oldIsExternalPrec, &oldApproxFwdOpSrc, &oldSupportSolveUse );
 
  634   typedef Belos::LinearProblem<Scalar,MV_t,LO_t> LP_t;
 
  640     lp = 
rcp(
new LP_t());
 
  647   lp->setOperator(fwdOp);
 
  658       !( left.
get() || right.
get() || unspecified.
get() ), std::logic_error
 
  659       ,
"Error, at least one preconditoner linear operator objects must be set!" 
  662       if (paramList_->get<
bool>(LeftPreconditionerIfUnspecified_name, 
false))
 
  663         lp->setLeftPrec(unspecified);
 
  665         lp->setRightPrec(unspecified);
 
  668       lp->setLeftPrec(left);
 
  671       lp->setRightPrec(right);
 
  677         ,
"Error, we can not currently handle split preconditioners!" 
  682     set_extra_data<RCP<PreconditionerBase<Scalar> > >(myPrec,
"Belos::InternalPrec",
 
  685   else if(prec.
get()) {
 
  686     set_extra_data<RCP<const PreconditionerBase<Scalar> > >(prec,
"Belos::ExternalPrec",
 
  694   typedef Belos::SolverManager<Scalar,MV_t,LO_t> IterativeSolver_t;
 
  698   switch(solverType_) {
 
  699     case SOLVER_TYPE_BLOCK_GMRES: 
 
  702       if(paramList_.get()) {
 
  709         iterativeSolver = oldIterSolver;
 
  710         iterativeSolver->setProblem( lp );
 
  711         iterativeSolver->setParameters( solverPL );
 
  714         iterativeSolver = 
rcp(
new Belos::BlockGmresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
 
  718     case SOLVER_TYPE_PSEUDO_BLOCK_GMRES:
 
  721       if(paramList_.get()) {
 
  730         iterativeSolver = oldIterSolver;
 
  731         iterativeSolver->setProblem( lp );
 
  732         iterativeSolver->setParameters( solverPL );
 
  735         iterativeSolver = 
rcp(
new Belos::PseudoBlockGmresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
 
  739     case SOLVER_TYPE_BLOCK_CG:
 
  742       if(paramList_.get()) {
 
  749         iterativeSolver = oldIterSolver;
 
  750         iterativeSolver->setProblem( lp );
 
  751         iterativeSolver->setParameters( solverPL );
 
  754         iterativeSolver = 
rcp(
new Belos::BlockCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
 
  758     case SOLVER_TYPE_PSEUDO_BLOCK_CG:
 
  761       if(paramList_.get()) {
 
  770         iterativeSolver = oldIterSolver;
 
  771         iterativeSolver->setProblem( lp );
 
  772         iterativeSolver->setParameters( solverPL );
 
  775         iterativeSolver = 
rcp(
new Belos::PseudoBlockCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
 
  779     case SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG:
 
  782       if(paramList_.get()) {
 
  791         iterativeSolver = oldIterSolver;
 
  792         iterativeSolver->setProblem( lp );
 
  793         iterativeSolver->setParameters( solverPL );
 
  796         iterativeSolver = 
rcp(
new Belos::PseudoBlockStochasticCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
 
  800     case SOLVER_TYPE_GCRODR:
 
  803       if(paramList_.get()) {
 
  810         iterativeSolver = oldIterSolver;
 
  811         iterativeSolver->setProblem( lp );
 
  812         iterativeSolver->setParameters( solverPL );
 
  815         iterativeSolver = 
rcp(
new Belos::GCRODRSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
 
  819     case SOLVER_TYPE_RCG:
 
  822       if(paramList_.get()) {
 
  829         iterativeSolver = oldIterSolver;
 
  830         iterativeSolver->setProblem( lp );
 
  831         iterativeSolver->setParameters( solverPL );
 
  834         iterativeSolver = 
rcp(
new Belos::RCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
 
  838     case SOLVER_TYPE_MINRES:
 
  841       if(paramList_.get()) {
 
  848         iterativeSolver = oldIterSolver;
 
  849         iterativeSolver->setProblem( lp );
 
  850         iterativeSolver->setParameters( solverPL );
 
  853         iterativeSolver = 
rcp(
new Belos::MinresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
 
  857     case SOLVER_TYPE_TFQMR:
 
  860       if(paramList_.get()) {
 
  867         iterativeSolver = oldIterSolver;
 
  868         iterativeSolver->setProblem( lp );
 
  869         iterativeSolver->setParameters( solverPL );
 
  872         iterativeSolver = 
rcp(
new Belos::TFQMRSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
 
  888     lp, solverPL, iterativeSolver,
 
  889     fwdOpSrc, prec, myPrec.
get()==NULL, approxFwdOpSrc,
 
  890     supportSolveUse, convergenceTestFrequency_
 
  892   belosOp->setOStream(out);
 
  893   belosOp->setVerbLevel(verbLevel);
 
  895   if(paramList_.get()) {
 
  897     paramList_->validateParameters(*this->getValidParameters(),1); 
 
  901     *out << 
"\nLeaving Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<
">::initializeOpImpl(...) ...\n";
 
  909 #endif // THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP 
void updateThisValidParamList()
 
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const 
 
void uninitializeOp(LinearOpWithSolveBase< Scalar > *Op, Teuchos::RCP< const LinearOpSourceBase< Scalar > > *fwdOpSrc, Teuchos::RCP< const PreconditionerBase< Scalar > > *prec, Teuchos::RCP< const LinearOpSourceBase< Scalar > > *approxFwdOpSrc, ESupportSolveUse *supportSolveUse) const 
 
RCP< const LinearOpSourceBase< Scalar > > extract_approxFwdOpSrc()
 
static void addConverter(RCP< const ParameterEntryValidator > validator, RCP< ValidatorXMLConverter > converterToAdd)
 
void initialize(const RCP< Belos::LinearProblem< Scalar, MV_t, LO_t > > &lp, const RCP< Teuchos::ParameterList > &solverPL, const RCP< Belos::SolverManager< Scalar, MV_t, LO_t > > &iterativeSolver, const RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const RCP< const PreconditionerBase< Scalar > > &prec, const bool isExternalPrec, const RCP< const LinearOpSourceBase< Scalar > > &approxFwdOpSrc, const ESupportSolveUse &supportSolveUse, const int convergenceTestFrequency)
Initializes given precreated solver objects. 
 
Concrete LinearOpWithSolveBase subclass in terms of Belos. 
 
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
RCP< const PreconditionerBase< Scalar > > extract_prec()
 
void initializeAndReuseOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, LinearOpWithSolveBase< Scalar > *Op) const 
 
RCP< ParameterList > sublist(const RCP< ParameterList > ¶mList, const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
 
static Teuchos::RCP< const Teuchos::ParameterList > generateAndGetValidParameters()
 
bool isExternalPrec() const 
 
void unsetPreconditionerFactory(Teuchos::RCP< PreconditionerFactoryBase< Scalar > > *precFactory, std::string *precFactoryName)
 
Thyra specializations of MultiVecTraits and OperatorTraits. 
 
void setParameterList(Teuchos::RCP< Teuchos::ParameterList > const ¶mList)
 
bool isCompatible(const LinearOpSourceBase< Scalar > &fwdOpSrc) const 
 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
 
bool acceptsPreconditionerFactory() const 
Returns true . 
 
T_To & dyn_cast(T_From &from)
 
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
 
void initializeOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const 
 
std::string description() const 
 
void validateParametersAndSetDefaults(ParameterList const &validParamList, int const depth=1000)
 
void setPreconditionerFactory(const Teuchos::RCP< PreconditionerFactoryBase< Scalar > > &precFactory, const std::string &precFactoryName)
 
TypeTo as(const TypeFrom &t)
 
ESupportSolveUse supportSolveUse() const 
 
bool supportsPreconditionerInputType(const EPreconditionerInputType precOpType) const 
 
Teuchos::RCP< const Teuchos::ParameterList > getParameterList() const 
 
RCP< const LinearOpSourceBase< Scalar > > extract_fwdOpSrc()
 
bool nonnull(const boost::shared_ptr< T > &p)
 
Teuchos::RCP< PreconditionerFactoryBase< Scalar > > getPreconditionerFactory() const 
 
void initializeApproxPreconditionedOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &approxFwdOpSrc, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const 
 
void initializePreconditionedOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const Teuchos::RCP< const PreconditionerBase< Scalar > > &prec, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const 
 
Teuchos::RCP< LinearOpWithSolveBase< Scalar > > createOp() const 
 
void uninitialize(RCP< Belos::LinearProblem< Scalar, MV_t, LO_t > > *lp=NULL, RCP< Teuchos::ParameterList > *solverPL=NULL, RCP< Belos::SolverManager< Scalar, MV_t, LO_t > > *iterativeSolver=NULL, RCP< const LinearOpSourceBase< Scalar > > *fwdOpSrc=NULL, RCP< const PreconditionerBase< Scalar > > *prec=NULL, bool *isExternalPrec=NULL, RCP< const LinearOpSourceBase< Scalar > > *approxFwdOpSrc=NULL, ESupportSolveUse *supportSolveUse=NULL)
Uninitializes and returns stored quantities. 
 
BelosLinearOpWithSolveFactory()
Construct without preconditioner factory. 
 
void initializeOpImpl(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &approxFwdOpSrc, const Teuchos::RCP< const PreconditionerBase< Scalar > > &prec, const bool reusePrec, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const 
 
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
 
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)