10 #ifndef THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP
11 #define THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP
15 #include "Thyra_BelosLinearOpWithSolve.hpp"
16 #include "Thyra_ScaledAdjointLinearOpBase.hpp"
18 #include "BelosBlockGmresSolMgr.hpp"
19 #include "BelosPseudoBlockGmresSolMgr.hpp"
20 #include "BelosBlockCGSolMgr.hpp"
21 #include "BelosPseudoBlockCGSolMgr.hpp"
22 #include "BelosPseudoBlockStochasticCGSolMgr.hpp"
23 #include "BelosGCRODRSolMgr.hpp"
24 #include "BelosRCGSolMgr.hpp"
25 #include "BelosMinresSolMgr.hpp"
26 #include "BelosTFQMRSolMgr.hpp"
27 #include "BelosBiCGStabSolMgr.hpp"
28 #include "BelosFixedPointSolMgr.hpp"
31 #if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
48 template<
class Scalar>
49 const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverType_name =
"Solver Type";
50 template<
class Scalar>
51 const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverType_default =
"Pseudo Block GMRES";
52 template<
class Scalar>
53 const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverTypes_name =
"Solver Types";
54 template<
class Scalar>
55 const std::string BelosLinearOpWithSolveFactory<Scalar>::BlockGMRES_name =
"Block GMRES";
56 template<
class Scalar>
57 const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockGMRES_name =
"Pseudo Block GMRES";
58 template<
class Scalar>
59 const std::string BelosLinearOpWithSolveFactory<Scalar>::BlockCG_name =
"Block CG";
60 template<
class Scalar>
61 const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockCG_name =
"Pseudo Block CG";
62 template<
class Scalar>
63 const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockStochasticCG_name =
"Pseudo Block Stochastic CG";
64 template<
class Scalar>
65 const std::string BelosLinearOpWithSolveFactory<Scalar>::GCRODR_name =
"GCRODR";
66 template<
class Scalar>
67 const std::string BelosLinearOpWithSolveFactory<Scalar>::RCG_name =
"RCG";
68 template<
class Scalar>
69 const std::string BelosLinearOpWithSolveFactory<Scalar>::MINRES_name =
"MINRES";
70 template<
class Scalar>
71 const std::string BelosLinearOpWithSolveFactory<Scalar>::TFQMR_name =
"TFQMR";
72 template<
class Scalar>
73 const std::string BelosLinearOpWithSolveFactory<Scalar>::BiCGStab_name =
"BiCGStab";
74 template<
class Scalar>
75 const std::string BelosLinearOpWithSolveFactory<Scalar>::FixedPoint_name =
"Fixed Point";
77 #if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
78 template<
class Scalar>
79 const std::string BelosLinearOpWithSolveFactory<Scalar>::TpetraGmres_name =
"TPETRA GMRES";
80 template<
class Scalar>
81 const std::string BelosLinearOpWithSolveFactory<Scalar>::TpetraGmresPipeline_name =
"TPETRA GMRES PIPELINE";
82 template<
class Scalar>
83 const std::string BelosLinearOpWithSolveFactory<Scalar>::TpetraGmresSingleReduce_name =
"TPETRA GMRES SINGLE REDUCE";
84 template<
class Scalar>
85 const std::string BelosLinearOpWithSolveFactory<Scalar>::TpetraGmresSstep_name =
"TPETRA GMRES S-STEP";
88 template<
class Scalar>
89 const std::string BelosLinearOpWithSolveFactory<Scalar>::ConvergenceTestFrequency_name =
"Convergence Test Frequency";
92 const std::string LeftPreconditionerIfUnspecified_name =
"Left Preconditioner If Unspecified";
98 template<
class Scalar>
100 :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES),
101 convergenceTestFrequency_(1)
107 template<
class Scalar>
109 const RCP<PreconditionerFactoryBase<Scalar> > &precFactory
111 :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES)
120 template<
class Scalar>
127 template<
class Scalar>
129 const RCP<PreconditionerFactoryBase<Scalar> > &precFactory,
130 const std::string &precFactoryName
135 precFactoryValidPL = precFactory->getValidParameters();
136 const std::string _precFactoryName =
137 ( precFactoryName !=
""
139 : ( precFactoryValidPL.get() ? precFactoryValidPL->name() :
"GENERIC PRECONDITIONER FACTORY" )
141 precFactory_ = precFactory;
142 precFactoryName_ = _precFactoryName;
143 updateThisValidParamList();
147 template<
class Scalar>
155 template<
class Scalar>
157 RCP<PreconditionerFactoryBase<Scalar> > *precFactory,
158 std::string *precFactoryName
161 if(precFactory) *precFactory = precFactory_;
162 if(precFactoryName) *precFactoryName = precFactoryName_;
164 precFactoryName_ =
"";
165 updateThisValidParamList();
169 template<
class Scalar>
171 const LinearOpSourceBase<Scalar> &fwdOpSrc
174 if(precFactory_.get())
175 return precFactory_->isCompatible(fwdOpSrc);
180 template<
class Scalar>
188 template<
class Scalar>
190 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
191 LinearOpWithSolveBase<Scalar> *Op,
192 const ESupportSolveUse supportSolveUse
196 initializeOpImpl(fwdOpSrc,null,null,
false,Op,supportSolveUse);
200 template<
class Scalar>
202 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
203 LinearOpWithSolveBase<Scalar> *Op
207 initializeOpImpl(fwdOpSrc,null,null,
true,Op,SUPPORT_SOLVE_UNSPECIFIED);
211 template<
class Scalar>
213 const EPreconditionerInputType precOpType
216 if(precFactory_.get())
218 return (precOpType==PRECONDITIONER_INPUT_TYPE_AS_OPERATOR);
222 template<
class Scalar>
224 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
225 const RCP<
const PreconditionerBase<Scalar> > &prec,
226 LinearOpWithSolveBase<Scalar> *Op,
227 const ESupportSolveUse supportSolveUse
231 initializeOpImpl(fwdOpSrc,null,prec,
false,Op,supportSolveUse);
235 template<
class Scalar>
237 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
238 const RCP<
const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
239 LinearOpWithSolveBase<Scalar> *Op,
240 const ESupportSolveUse supportSolveUse
244 initializeOpImpl(fwdOpSrc,approxFwdOpSrc,null,
false,Op,supportSolveUse);
248 template<
class Scalar>
250 LinearOpWithSolveBase<Scalar> *Op,
251 RCP<
const LinearOpSourceBase<Scalar> > *fwdOpSrc,
252 RCP<
const PreconditionerBase<Scalar> > *prec,
253 RCP<
const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
254 ESupportSolveUse *supportSolveUse
273 if(fwdOpSrc) *fwdOpSrc = _fwdOpSrc;
274 if(prec) *prec = _prec;
275 if(approxFwdOpSrc) *approxFwdOpSrc = _approxFwdOpSrc;
276 if(supportSolveUse) *supportSolveUse = _supportSolveUse;
283 template<
class Scalar>
290 paramList_ = paramList;
292 Teuchos::getIntegralValue<EBelosSolverType>(*paramList_, SolverType_name);
293 convergenceTestFrequency_ =
294 Teuchos::getParameter<int>(*paramList_, ConvergenceTestFrequency_name);
295 Teuchos::readVerboseObjectSublist(&*paramList_,
this);
299 template<
class Scalar>
307 template<
class Scalar>
317 template<
class Scalar>
325 template<
class Scalar>
329 return thisValidParamList_;
336 template<
class Scalar>
339 std::ostringstream oss;
343 oss << precFactory_->description();
353 template<
class Scalar>
358 using Teuchos::tuple;
359 using Teuchos::setStringToIntegralParameter;
365 EBelosSolverType> >::getDummyObject());
367 typedef MultiVectorBase<Scalar> MV_t;
368 typedef LinearOpBase<Scalar> LO_t;
370 if(validParamList.
get()==NULL) {
372 setStringToIntegralParameter<EBelosSolverType>(
373 SolverType_name, SolverType_default,
374 "Type of linear solver algorithm to use.",
377 "Pseudo Block GMRES",
380 "Pseudo Block Stochastic CG",
387 #if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
389 "TPETRA GMRES PIPELINE",
390 "TPETRA GMRES SINGLE REDUCE",
391 "TPETRA GMRES S-STEP"
395 "Block GMRES solver for nonsymmetric linear systems. It can also solve "
396 "single right-hand side systems, and can also perform Flexible GMRES "
397 "(where the preconditioner may change at every iteration, for example "
398 "for inner-outer iterations), by setting options in the \"Block GMRES\" "
401 "GMRES solver for nonsymmetric linear systems, that performs single "
402 "right-hand side solves on multiple right-hand sides at once. It "
403 "exploits operator multivector multiplication in order to amortize "
404 "global communication costs. Individual linear systems are deflated "
405 "out as they are solved.",
407 "Block CG solver for symmetric (Hermitian in complex arithmetic) "
408 "positive definite linear systems. It can also solve single "
409 "right-hand-side systems.",
411 "CG solver that performs single right-hand side CG on multiple right-hand "
412 "sides at once. It exploits operator multivector multiplication in order "
413 "to amortize global communication costs. Individual linear systems are "
414 "deflated out as they are solved.",
416 "stochastic CG solver that performs single right-hand side CG on multiple right-hand "
417 "sides at once. It exploits operator multivector multiplication in order "
418 "to amortize global communication costs. Individual linear systems are "
419 "deflated out as they are solved. [EXPERIMENTAL]",
421 "Variant of GMRES that performs subspace recycling to accelerate "
422 "convergence for sequences of solves with related linear systems. "
423 "Individual linear systems are deflated out as they are solved. "
424 "The current implementation only supports real-valued Scalar types.",
426 "CG solver for symmetric (Hermitian in complex arithmetic) positive "
427 "definite linear systems, that performs subspace recycling to "
428 "accelerate convergence for sequences of related linear systems.",
430 "MINRES solver for symmetric indefinite linear systems. It performs "
431 "single-right-hand-side solves on multiple right-hand sides sequentially.",
433 "TFQMR (Transpose-Free QMR) solver for nonsymmetric linear systems.",
435 "BiCGStab solver for nonsymmetric linear systems.",
437 "Fixed point iteration"
439 #
if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
440 ,
"Native Tpetra implementation of GMRES",
442 "Native Tpetra implementation of pipeline GMRES",
444 "Native Tpetra implementation of single-reduce GMRES",
446 "Native Tpetra implementation of s-step GMRES"
449 tuple<EBelosSolverType>(
450 SOLVER_TYPE_BLOCK_GMRES,
451 SOLVER_TYPE_PSEUDO_BLOCK_GMRES,
452 SOLVER_TYPE_BLOCK_CG,
453 SOLVER_TYPE_PSEUDO_BLOCK_CG,
454 SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG,
459 SOLVER_TYPE_BICGSTAB,
460 SOLVER_TYPE_FIXEDPOINT
461 #
if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
462 ,SOLVER_TYPE_TPETRA_GMRES,
463 SOLVER_TYPE_TPETRA_GMRES_PIPELINE,
464 SOLVER_TYPE_TPETRA_GMRES_SINGLE_REDUCE,
465 SOLVER_TYPE_TPETRA_GMRES_SSTEP
470 validParamList->
set(ConvergenceTestFrequency_name, as<int>(1),
471 "Number of linear solver iterations to skip between applying"
472 " user-defined convergence test.");
474 LeftPreconditionerIfUnspecified_name,
false,
475 "If the preconditioner does not specify if it is left or right, and this\n"
476 "option is set to true, put the preconditioner on the left side.\n"
477 "Historically, preconditioning is on the right. Some solvers may not\n"
478 "support left preconditioning.");
480 &solverTypesSL = validParamList->
sublist(SolverTypes_name);
482 const bool lapackSupportsScalar = Belos::Details::LapackSupportsScalar<Scalar>::value;
486 Belos::BlockGmresSolMgr<Scalar,MV_t,LO_t> mgr;
487 solverTypesSL.
sublist(BlockGMRES_name).setParameters(
488 *mgr.getValidParameters()
492 Belos::PseudoBlockGmresSolMgr<Scalar,MV_t,LO_t> mgr;
493 solverTypesSL.
sublist(PseudoBlockGMRES_name).setParameters(
494 *mgr.getValidParameters()
497 if (lapackSupportsScalar) {
498 Belos::BlockCGSolMgr<Scalar,MV_t,LO_t> mgr;
499 solverTypesSL.
sublist(BlockCG_name).setParameters(
500 *mgr.getValidParameters()
503 if (lapackSupportsScalar) {
504 Belos::PseudoBlockCGSolMgr<Scalar,MV_t,LO_t> mgr;
505 solverTypesSL.
sublist(PseudoBlockCG_name).setParameters(
506 *mgr.getValidParameters()
510 Belos::PseudoBlockStochasticCGSolMgr<Scalar,MV_t,LO_t> mgr;
511 solverTypesSL.
sublist(PseudoBlockStochasticCG_name).setParameters(
512 *mgr.getValidParameters()
515 if (lapackSupportsScalar) {
516 Belos::GCRODRSolMgr<Scalar,MV_t,LO_t> mgr;
517 solverTypesSL.
sublist(GCRODR_name).setParameters(
518 *mgr.getValidParameters()
521 if (lapackSupportsScalar && scalarIsReal) {
522 Belos::RCGSolMgr<Scalar,MV_t,LO_t> mgr;
523 solverTypesSL.
sublist(RCG_name).setParameters(
524 *mgr.getValidParameters()
528 Belos::MinresSolMgr<Scalar,MV_t,LO_t> mgr;
529 solverTypesSL.
sublist(MINRES_name).setParameters(
530 *mgr.getValidParameters()
534 Belos::TFQMRSolMgr<Scalar,MV_t,LO_t> mgr;
535 solverTypesSL.
sublist(TFQMR_name).setParameters(
536 *mgr.getValidParameters()
540 Belos::BiCGStabSolMgr<Scalar,MV_t,LO_t> mgr;
541 solverTypesSL.
sublist(BiCGStab_name).setParameters(
542 *mgr.getValidParameters()
546 Belos::FixedPointSolMgr<Scalar,MV_t,LO_t> mgr;
547 solverTypesSL.
sublist(FixedPoint_name).setParameters(
548 *mgr.getValidParameters()
551 #if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
554 solverTypesSL.
sublist(TpetraGmres_name).setParameters(
560 solverTypesSL.
sublist(TpetraGmresPipeline_name).setParameters(
566 solverTypesSL.
sublist(TpetraGmresSingleReduce_name).setParameters(
572 solverTypesSL.
sublist(TpetraGmresSstep_name).setParameters(
578 return validParamList;
582 template<
class Scalar>
588 Teuchos::setupVerboseObjectSublist(&*thisValidParamList_);
592 template<
class Scalar>
594 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
595 const RCP<
const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
596 const RCP<
const PreconditionerBase<Scalar> > &prec_in,
597 const bool reusePrec,
598 LinearOpWithSolveBase<Scalar> *Op,
599 const ESupportSolveUse supportSolveUse
604 using Teuchos::set_extra_data;
606 typedef MultiVectorBase<Scalar> MV_t;
607 typedef LinearOpBase<Scalar> LO_t;
613 *out <<
"\nEntering Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<
">::initializeOpImpl(...) ...\n";
624 fwdOp = fwdOpSrc->getOp(),
625 approxFwdOp = ( approxFwdOpSrc.get() ? approxFwdOpSrc->getOp() :
Teuchos::null );
646 if(precFactory_.get()) {
649 ? Teuchos::rcp_const_cast<PreconditionerBase<Scalar> >(belosOp->
extract_prec())
652 bool hasExistingPrec =
false;
654 hasExistingPrec =
true;
659 hasExistingPrec =
false;
660 myPrec = precFactory_->createPrec();
662 if( hasExistingPrec && reusePrec ) {
667 if(approxFwdOp.get())
668 precFactory_->initializePrec(approxFwdOpSrc,&*myPrec);
670 precFactory_->initializePrec(fwdOpSrc,&*myPrec);
680 bool oldIsExternalPrec =
false;
685 ESupportSolveUse oldSupportSolveUse = SUPPORT_SOLVE_UNSPECIFIED;
687 belosOp->
uninitialize( &oldLP, NULL, &oldIterSolver, &oldFwdOpSrc,
688 NULL, &oldIsExternalPrec, &oldApproxFwdOpSrc, &oldSupportSolveUse );
695 typedef Belos::LinearProblem<Scalar,MV_t,LO_t> LP_t;
701 lp =
rcp(
new LP_t());
708 lp->setOperator(fwdOp);
719 !( left.
get() || right.
get() || unspecified.
get() ), std::logic_error
720 ,
"Error, at least one preconditoner linear operator objects must be set!"
723 if (paramList_->get<
bool>(LeftPreconditionerIfUnspecified_name,
false))
724 lp->setLeftPrec(unspecified);
726 lp->setRightPrec(unspecified);
729 lp->setLeftPrec(left);
732 lp->setRightPrec(right);
738 ,
"Error, we can not currently handle split preconditioners!"
743 set_extra_data<RCP<PreconditionerBase<Scalar> > >(myPrec,
"Belos::InternalPrec",
746 else if(prec.
get()) {
747 set_extra_data<RCP<const PreconditionerBase<Scalar> > >(prec,
"Belos::ExternalPrec",
755 typedef Belos::SolverManager<Scalar,MV_t,LO_t> IterativeSolver_t;
759 switch(solverType_) {
760 case SOLVER_TYPE_BLOCK_GMRES:
763 if(paramList_.get()) {
770 iterativeSolver = oldIterSolver;
771 iterativeSolver->setProblem( lp );
772 iterativeSolver->setParameters( solverPL );
775 iterativeSolver =
rcp(
new Belos::BlockGmresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
779 case SOLVER_TYPE_PSEUDO_BLOCK_GMRES:
782 if(paramList_.get()) {
791 iterativeSolver = oldIterSolver;
792 iterativeSolver->setProblem( lp );
793 iterativeSolver->setParameters( solverPL );
796 iterativeSolver =
rcp(
new Belos::PseudoBlockGmresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
800 case SOLVER_TYPE_BLOCK_CG:
803 if(paramList_.get()) {
810 iterativeSolver = oldIterSolver;
811 iterativeSolver->setProblem( lp );
812 iterativeSolver->setParameters( solverPL );
815 iterativeSolver =
rcp(
new Belos::BlockCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
819 case SOLVER_TYPE_PSEUDO_BLOCK_CG:
822 if(paramList_.get()) {
831 iterativeSolver = oldIterSolver;
832 iterativeSolver->setProblem( lp );
833 iterativeSolver->setParameters( solverPL );
836 iterativeSolver =
rcp(
new Belos::PseudoBlockCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
840 case SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG:
843 if(paramList_.get()) {
852 iterativeSolver = oldIterSolver;
853 iterativeSolver->setProblem( lp );
854 iterativeSolver->setParameters( solverPL );
857 iterativeSolver =
rcp(
new Belos::PseudoBlockStochasticCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
861 case SOLVER_TYPE_GCRODR:
864 if(paramList_.get()) {
871 iterativeSolver = oldIterSolver;
872 iterativeSolver->setProblem( lp );
873 iterativeSolver->setParameters( solverPL );
876 iterativeSolver =
rcp(
new Belos::GCRODRSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
880 case SOLVER_TYPE_RCG:
883 if(paramList_.get()) {
890 iterativeSolver = oldIterSolver;
891 iterativeSolver->setProblem( lp );
892 iterativeSolver->setParameters( solverPL );
895 iterativeSolver =
rcp(
new Belos::RCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
899 case SOLVER_TYPE_MINRES:
902 if(paramList_.get()) {
909 iterativeSolver = oldIterSolver;
910 iterativeSolver->setProblem( lp );
911 iterativeSolver->setParameters( solverPL );
914 iterativeSolver =
rcp(
new Belos::MinresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
918 case SOLVER_TYPE_TFQMR:
921 if(paramList_.get()) {
928 iterativeSolver = oldIterSolver;
929 iterativeSolver->setProblem( lp );
930 iterativeSolver->setParameters( solverPL );
933 iterativeSolver =
rcp(
new Belos::TFQMRSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
937 case SOLVER_TYPE_BICGSTAB:
940 if(paramList_.get()) {
947 iterativeSolver = oldIterSolver;
948 iterativeSolver->setProblem( lp );
949 iterativeSolver->setParameters( solverPL );
952 iterativeSolver =
rcp(
new Belos::BiCGStabSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
956 case SOLVER_TYPE_FIXEDPOINT:
959 if(paramList_.get()) {
966 iterativeSolver = oldIterSolver;
967 iterativeSolver->setProblem( lp );
968 iterativeSolver->setParameters( solverPL );
971 iterativeSolver =
rcp(
new Belos::FixedPointSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
975 #if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
976 case SOLVER_TYPE_TPETRA_GMRES:
979 if(paramList_.get()) {
987 iterativeSolver->setProblem( lp );
988 iterativeSolver->setParameters( solverPL );
991 case SOLVER_TYPE_TPETRA_GMRES_PIPELINE:
994 if(paramList_.get()) {
997 solverPL =
Teuchos::rcp( &tpetraGmresPipelinePL,
false );
1002 iterativeSolver->setProblem( lp );
1003 iterativeSolver->setParameters( solverPL );
1006 case SOLVER_TYPE_TPETRA_GMRES_SINGLE_REDUCE:
1009 if(paramList_.get()) {
1012 solverPL =
Teuchos::rcp( &tpetraGmresSingleReducePL,
false );
1017 iterativeSolver->setProblem( lp );
1018 iterativeSolver->setParameters( solverPL );
1021 case SOLVER_TYPE_TPETRA_GMRES_SSTEP:
1024 if(paramList_.get()) {
1032 iterativeSolver->setProblem( lp );
1033 iterativeSolver->setParameters( solverPL );
1048 lp, solverPL, iterativeSolver,
1049 fwdOpSrc, prec, myPrec.
get()==NULL, approxFwdOpSrc,
1050 supportSolveUse, convergenceTestFrequency_
1052 belosOp->setOStream(out);
1053 belosOp->setVerbLevel(verbLevel);
1054 #ifdef TEUCHOS_DEBUG
1055 if(paramList_.get()) {
1057 paramList_->validateParameters(*this->getValidParameters(),1);
1061 *out <<
"\nLeaving Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<
">::initializeOpImpl(...) ...\n";
1069 #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.
#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="")
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
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::RCP< const Teuchos::ParameterList > getValidParameters() const override
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()
virtual std::string description() const
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)